GRASS Programmer's Manual  6.5.svn(2014)-r66266
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
N_pde.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * MODULE: Grass PDE Numerical Library
5 * AUTHOR(S): Soeren Gebbert, Berlin (GER) Dec 2006
6 * soerengebbert <at> gmx <dot> de
7 *
8 * PURPOSE: This file contains definitions of variables and data types
9 *
10 * COPYRIGHT: (C) 2000 by the GRASS Development Team
11 *
12 * This program is free software under the GNU General Public
13 * License (>=v2). Read the file COPYING that comes with GRASS
14 * for details.
15 *
16 *****************************************************************************/
17 
18 #include <grass/gis.h>
19 #include <grass/G3d.h>
20 #include <grass/glocale.h>
21 #include <grass/gmath.h>
22 
23 #ifndef _N_PDE_H_
24 #define _N_PDE_H_
25 
26 /*solver names */
27 #define N_SOLVER_DIRECT_GAUSS "gauss"
28 #define N_SOLVER_DIRECT_LU "lu"
29 #define N_SOLVER_DIRECT_CHOLESKY "cholesky"
30 #define N_SOLVER_ITERATIVE_JACOBI "jacobi"
31 #define N_SOLVER_ITERATIVE_SOR "sor"
32 #define N_SOLVER_ITERATIVE_CG "cg"
33 #define N_SOLVER_ITERATIVE_PCG "pcg"
34 #define N_SOLVER_ITERATIVE_BICGSTAB "bicgstab"
35 
36 /*preconditioner */
37 #define N_DIAGONAL_PRECONDITION 1
38 #define N_ROWSCALE_ABSSUMNORM_PRECONDITION 2
39 #define N_ROWSCALE_EUKLIDNORM_PRECONDITION 3
40 #define N_ROWSCALE_MAXNORM_PRECONDITION 4
41 
42 #define N_NORMAL_LES 0
43 #define N_SPARSE_LES 1
44 
45 #define N_CELL_INACTIVE 0
46 #define N_CELL_ACTIVE 1
47 #define N_CELL_DIRICHLET 2
48 #define N_CELL_TRANSMISSION 3
49 
52 #define N_MAX_CELL_STATE 20
53 
54 #define N_5_POINT_STAR 0
55 #define N_7_POINT_STAR 1
56 #define N_9_POINT_STAR 2
57 #define N_27_POINT_STAR 3
58 
59 #define N_MAXIMUM_NORM 0
60 #define N_EUKLID_NORM 1
61 
62 #define N_ARRAY_SUM 0 /* summ two arrays */
63 #define N_ARRAY_DIF 1 /* calc the difference between two arrays */
64 #define N_ARRAY_MUL 2 /* multiply two arrays */
65 #define N_ARRAY_DIV 3 /* array division, if div with 0 the NULL value is set */
66 
67 #define N_UPWIND_FULL 0 /*full upwinding stabilization */
68 #define N_UPWIND_EXP 1 /*exponential upwinding stabilization */
69 #define N_UPWIND_WEIGHT 2 /*weighted upwinding stabilization */
70 
71 
72 
73 /* *************************************************************** */
74 /* *************** LINEARE EQUATION SYSTEM PART ****************** */
75 /* *************************************************************** */
76 
86 typedef struct
87 {
88  double *x; /*the value vector */
89  double *b; /*the right side of Ax = b */
90  double **A; /*the normal quadratic matrix */
91  G_math_spvector **Asp; /*the sparse matrix */
92  int rows; /*number of rows */
93  int cols; /*number of cols */
94  int quad; /*is the matrix quadratic (1-quadratic, 0 not) */
95  int type; /*the type of the les, normal == 0, sparse == 1 */
96 } N_les;
97 
98 extern N_les *N_alloc_les_param(int cols, int rows, int type, int param);
99 
100 extern N_les *N_alloc_les(int rows, int type);
101 
102 extern N_les *N_alloc_les_A(int rows, int type);
103 
104 extern N_les *N_alloc_les_Ax(int rows, int type);
105 
106 extern N_les *N_alloc_les_Ax_b(int rows, int type);
107 
108 extern N_les *N_alloc_nquad_les(int cols, int rows, int type);
109 
110 extern N_les *N_alloc_nquad_les_A(int cols, int rows, int type);
111 
112 extern N_les *N_alloc_nquad_les_Ax(int cols, int rows, int type);
113 
114 extern N_les *N_alloc_nquad_les_Ax_b(int cols, int rows, int type);
115 
116 extern void N_print_les(N_les * les);
117 
118 extern void N_free_les(N_les * les);
119 
120 /* *************************************************************** */
121 /* *************** GEOMETRY INFORMATION ************************** */
122 /* *************************************************************** */
123 
127 typedef struct
128 {
129  int planimetric; /*If the projection is not planimetric (0), the array calculation is different for each row */
130  double *area; /* the vector of area values for non-planimetric projection for each row */
131  int dim; /* 2 or 3 */
132 
133  double dx;
134  double dy;
135  double dz;
136 
137  double Az;
138 
139  int depths;
140  int rows;
141  int cols;
142 
143 } N_geom_data;
144 
145 extern N_geom_data *N_alloc_geom_data(void);
146 
147 extern void N_free_geom_data(N_geom_data * geodata);
148 
149 extern N_geom_data *N_init_geom_data_3d(G3D_Region * region3d,
150  N_geom_data * geodata);
151 extern N_geom_data *N_init_geom_data_2d(struct Cell_head *region,
152  N_geom_data * geodata);
153 extern double N_get_geom_data_area_of_cell(N_geom_data * geom, int row);
154 
155 /* *************************************************************** */
156 /* *************** READING RASTER AND VOLUME DATA **************** */
157 /* *************************************************************** */
158 
159 typedef struct
160 {
161  int type; /* which raster type CELL_TYPE, FCELL_TYPE, DCELL_TYPE */
162  int rows, cols;
163  int rows_intern, cols_intern;
164  int offset; /*number of cols/rows offset at each boundary */
165  CELL *cell_array; /*The data is stored in an one dimensional array internally */
166  FCELL *fcell_array; /*The data is stored in an one dimensional array internally */
167  DCELL *dcell_array; /*The data is stored in an one dimensional array internally */
168 } N_array_2d;
169 
170 extern N_array_2d *N_alloc_array_2d(int cols, int rows, int offset, int type);
171 
172 extern void N_free_array_2d(N_array_2d * data_array);
173 
174 extern int N_get_array_2d_type(N_array_2d * array2d);
175 
176 extern void N_get_array_2d_value(N_array_2d * array2d, int col, int row,
177  void *value);
178 extern CELL N_get_array_2d_c_value(N_array_2d * array2d, int col, int row);
179 
180 extern FCELL N_get_array_2d_f_value(N_array_2d * array2d, int col, int row);
181 
182 extern DCELL N_get_array_2d_d_value(N_array_2d * array2d, int col, int row);
183 
184 extern void N_put_array_2d_value(N_array_2d * array2d, int col, int row,
185  char *value);
186 extern void N_put_array_2d_c_value(N_array_2d * array2d, int col, int row,
187  CELL value);
188 extern void N_put_array_2d_f_value(N_array_2d * array2d, int col, int row,
189  FCELL value);
190 extern void N_put_array_2d_d_value(N_array_2d * array2d, int col, int row,
191  DCELL value);
192 extern int N_is_array_2d_value_null(N_array_2d * array2d, int col, int row);
193 
194 extern void N_put_array_2d_value_null(N_array_2d * array2d, int col, int row);
195 
196 extern void N_print_array_2d(N_array_2d * data);
197 
198 extern void N_print_array_2d_info(N_array_2d * data);
199 
200 extern void N_copy_array_2d(N_array_2d * source, N_array_2d * target);
201 
202 extern double N_norm_array_2d(N_array_2d * array1, N_array_2d * array2,
203  int type);
204 extern N_array_2d *N_math_array_2d(N_array_2d * array1, N_array_2d * array2,
205  N_array_2d * result, int type);
207 
208 extern N_array_2d *N_read_rast_to_array_2d(char *name, N_array_2d * array);
209 
210 extern void N_write_array_2d_to_rast(N_array_2d * array, char *name);
211 
212 extern void N_calc_array_2d_stats(N_array_2d * a, double *min, double *max,
213  double *sum, int *nonzero, int withoffset);
214 
215 typedef struct
216 {
217  int type; /* which raster type FCELL_TYPE, DCELL_TYPE */
218  int rows, cols, depths;
219  int rows_intern, cols_intern, depths_intern;
220  int offset; /*number of cols/rows/depths offset at each boundary */
221  float *fcell_array; /*The data is stored in an one dimensional array internally */
222  double *dcell_array; /*The data is stored in an one dimensional array internally */
223 } N_array_3d;
224 
225 extern N_array_3d *N_alloc_array_3d(int cols, int rows, int depths,
226  int offset, int type);
227 extern void N_free_array_3d(N_array_3d * data_array);
228 
229 extern int N_get_array_3d_type(N_array_3d * array3d);
230 
231 extern void N_get_array_3d_value(N_array_3d * array3d, int col, int row,
232  int depth, void *value);
233 extern float N_get_array_3d_f_value(N_array_3d * array3d, int col, int row,
234  int depth);
235 extern double N_get_array_3d_d_value(N_array_3d * array3d, int col, int row,
236  int depth);
237 extern void N_put_array_3d_value(N_array_3d * array3d, int col, int row,
238  int depth, char *value);
239 extern void N_put_array_3d_f_value(N_array_3d * array3d, int col, int row,
240  int depth, float value);
241 extern void N_put_array_3d_d_value(N_array_3d * array3d, int col, int row,
242  int depth, double value);
243 extern int N_is_array_3d_value_null(N_array_3d * array3d, int col, int row,
244  int depth);
245 extern void N_put_array_3d_value_null(N_array_3d * array3d, int col, int row,
246  int depth);
247 extern void N_print_array_3d(N_array_3d * data);
248 
249 extern void N_print_array_3d_info(N_array_3d * data);
250 
251 extern void N_copy_array_3d(N_array_3d * source, N_array_3d * target);
252 
253 extern double N_norm_array_3d(N_array_3d * array1, N_array_3d * array2,
254  int type);
255 extern N_array_3d *N_math_array_3d(N_array_3d * array1, N_array_3d * array2,
256  N_array_3d * result, int type);
258 
259 extern N_array_3d *N_read_rast3d_to_array_3d(char *name, N_array_3d * array,
260  int mask);
261 extern void N_write_array_3d_to_rast3d(N_array_3d * array, char *name,
262  int mask);
263 extern void N_calc_array_3d_stats(N_array_3d * a, double *min, double *max,
264  double *sum, int *nonzero, int withoffset);
265 
266 /* *************************************************************** */
267 /* *************** MATRIX ASSEMBLING METHODS ********************* */
268 /* *************************************************************** */
342 typedef struct
343 {
344  int type;
345  int count;
346  double C, W, E, N, S, NE, NW, SE, SW, V;
347  /*top part */
348  double T, W_T, E_T, N_T, S_T, NE_T, NW_T, SE_T, SW_T;
349  /*bottom part */
350  double B, W_B, E_B, N_B, S_B, NE_B, NW_B, SE_B, SW_B;
351 } N_data_star;
352 
356 typedef struct
357 {
358  N_data_star *(*callback) ();
360 
364 typedef struct
365 {
366  N_data_star *(*callback) ();
368 
369 
371  N_data_star * (*callback_func_3d) ());
373  N_data_star * (*callback_func_2d) ());
375 
377 
378 extern N_data_star *N_alloc_5star(void);
379 
380 extern N_data_star *N_alloc_7star(void);
381 
382 extern N_data_star *N_alloc_9star(void);
383 
384 extern N_data_star *N_alloc_27star(void);
385 
386 extern N_data_star *N_create_5star(double C, double W, double E, double N,
387  double S, double V);
388 extern N_data_star *N_create_7star(double C, double W, double E, double N,
389  double S, double T, double B, double V);
390 extern N_data_star *N_create_9star(double C, double W, double E, double N,
391  double S, double NW, double SW, double NE,
392  double SE, double V);
393 extern N_data_star *N_create_27star(double C, double W, double E, double N,
394  double S, double NW, double SW, double NE,
395  double SE, double T, double W_T,
396  double E_T, double N_T, double S_T,
397  double NW_T, double SW_T, double NE_T,
398  double SE_T, double B, double W_B,
399  double E_B, double N_B, double S_B,
400  double NW_B, double SW_B, double NE_B,
401  double SE_B, double V);
402 
403 extern N_data_star *N_callback_template_3d(void *data, N_geom_data * geom,
404  int col, int row, int depth);
405 extern N_data_star *N_callback_template_2d(void *data, N_geom_data * geom,
406  int col, int row);
407 extern N_les *N_assemble_les_3d(int les_type, N_geom_data * geom,
408  N_array_3d * status, N_array_3d * start_val,
409  void *data, N_les_callback_3d * callback);
410 extern N_les *N_assemble_les_3d_active(int les_type, N_geom_data * geom,
411  N_array_3d * status,
412  N_array_3d * start_val, void *data,
413  N_les_callback_3d * callback);
414 extern N_les *N_assemble_les_3d_dirichlet(int les_type, N_geom_data * geom,
415  N_array_3d * status,
416  N_array_3d * start_val, void *data,
417  N_les_callback_3d * callback);
418 extern N_les *N_assemble_les_3d_param(int les_type, N_geom_data * geom,
419  N_array_3d * status,
420  N_array_3d * start_val, void *data,
421  N_les_callback_3d * callback,
422  int cell_type);
423 extern N_les *N_assemble_les_2d(int les_type, N_geom_data * geom,
424  N_array_2d * status, N_array_2d * start_val,
425  void *data, N_les_callback_2d * callback);
426 extern N_les *N_assemble_les_2d_active(int les_type, N_geom_data * geom,
427  N_array_2d * status,
428  N_array_2d * start_val, void *data,
429  N_les_callback_2d * callback);
430 extern N_les *N_assemble_les_2d_dirichlet(int les_type, N_geom_data * geom,
431  N_array_2d * status,
432  N_array_2d * start_val, void *data,
433  N_les_callback_2d * callback);
434 extern N_les *N_assemble_les_2d_param(int les_type, N_geom_data * geom,
435  N_array_2d * status,
436  N_array_2d * start_val, void *data,
437  N_les_callback_2d * callback,
438  int cell_Type);
439 
440 extern int N_les_pivot_create(N_les * les);
441 
443  N_array_2d * status, N_array_2d * start_val);
445  N_array_3d * status, N_array_3d * start_val);
446 
447 /* *************************************************************** */
448 /* *************** GPDE STANDARD OPTIONS ************************* */
449 /* *************************************************************** */
450 
453 typedef enum
454 {
461 } N_STD_OPT;
462 
463 extern struct Option *N_define_standard_option(int opt);
464 
465 /* *************************************************************** */
466 /* *************** GPDE MATHEMATICAL TOOLS *********************** */
467 /* *************************************************************** */
468 
469 extern double N_calc_arith_mean(double a, double b);
470 
471 extern double N_calc_arith_mean_n(double *a, int size);
472 
473 extern double N_calc_geom_mean(double a, double b);
474 
475 extern double N_calc_geom_mean_n(double *a, int size);
476 
477 extern double N_calc_harmonic_mean(double a, double b);
478 
479 extern double N_calc_harmonic_mean_n(double *a, int size);
480 
481 extern double N_calc_quad_mean(double a, double b);
482 
483 extern double N_calc_quad_mean_n(double *a, int size);
484 
485 /* *************************************************************** */
486 /* *************** UPWIND STABILIZATION ALGORITHMS *************** */
487 /* *************************************************************** */
488 
489 extern double N_full_upwinding(double sprod, double distance, double D);
490 
491 extern double N_exp_upwinding(double sprod, double distance, double D);
492 
493 
494 /* *************************************************************** */
495 /* *************** METHODS FOR GRADIENT CALCULATION ************** */
496 /* *************************************************************** */
525 typedef struct
526 {
527 
528  double NC, SC, WC, EC;
529 
530 } N_gradient_2d;
531 
533 typedef struct
534 {
535 
536  double NC, SC, WC, EC, TC, BC;
537 
538 } N_gradient_3d;
539 
540 
585 typedef struct
586 {
587 
588  double NWN, NEN, WC, EC, SWS, SES;
589 
591 
593 typedef struct
594 {
595 
596  double NWW, NEE, NC, SC, SWW, SEE;
597 
599 
601 typedef struct
602 {
603 
604  double NWZ, NZ, NEZ, WZ, CZ, EZ, SWZ, SZ, SEZ;
605 
607 
609 typedef struct
610 {
611 
614 
616 
617 
619 typedef struct
620 {
621 
622  N_gradient_neighbours_x *xt; /*top values */
623  N_gradient_neighbours_x *xc; /*center values */
624  N_gradient_neighbours_x *xb; /*bottom values */
625 
626  N_gradient_neighbours_y *yt; /*top values */
627  N_gradient_neighbours_y *yc; /*center values */
628  N_gradient_neighbours_y *yb; /*bottom values */
629 
630  N_gradient_neighbours_z *zt; /*top-center values */
631  N_gradient_neighbours_z *zb; /*bottom-center values */
632 
634 
635 
637 typedef struct
638 {
639 
642  int cols, rows;
643  double min, max, mean, sum;
644  int nonull;
645 
647 
649 typedef struct
650 {
651 
655  int cols, rows, depths;
656  double min, max, mean, sum;
657  int nonull;
658 
660 
661 
662 extern N_gradient_2d *N_alloc_gradient_2d(void);
663 
664 extern void N_free_gradient_2d(N_gradient_2d * grad);
665 
666 extern N_gradient_2d *N_create_gradient_2d(double NC, double SC, double WC,
667  double EC);
668 extern int N_copy_gradient_2d(N_gradient_2d * source, N_gradient_2d * target);
669 
671  N_gradient_2d * gradient, int col,
672  int row);
673 
674 extern N_gradient_3d *N_alloc_gradient_3d(void);
675 
676 extern void N_free_gradient_3d(N_gradient_3d * grad);
677 
678 extern N_gradient_3d *N_create_gradient_3d(double NC, double SC, double WC,
679  double EC, double TC, double BC);
680 extern int N_copy_gradient_3d(N_gradient_3d * source, N_gradient_3d * target);
681 
683  N_gradient_3d * gradient, int col,
684  int row, int depth);
685 
687 
689 
691  double NEN,
692  double WC,
693  double EC,
694  double SWS,
695  double SES);
697  N_gradient_neighbours_x * target);
698 
700 
702 
704  double NEE,
705  double NC,
706  double SC,
707  double SWW,
708  double SEE);
710  N_gradient_neighbours_y * target);
711 
713 
715 
717  double NZ,
718  double NEZ,
719  double WZ,
720  double CZ,
721  double EZ,
722  double SWZ,
723  double SZ,
724  double SEZ);
726  N_gradient_neighbours_z * target);
727 
729 
731 
736  N_gradient_neighbours_2d * target);
739  N_gradient_neighbours_2d * gradient,
740  int col, int row);
741 
742 
744 
746 
757  N_gradient_neighbours_3d * target);
758 
760 
762 
763 
764 extern N_gradient_field_2d *N_alloc_gradient_field_2d(int cols, int rows);
765 
766 extern void N_free_gradient_field_2d(N_gradient_field_2d * field);
767 
769  N_gradient_field_2d * target);
771  N_array_2d * weight_x,
772  N_array_2d * weight_y,
773  N_geom_data * geom,
775  gradfield);
777  field, N_array_2d * x_comp,
778  N_array_2d * y_comp);
779 
781 
783 
785  int depths);
786 extern void N_free_gradient_field_3d(N_gradient_field_3d * field);
787 
789  N_gradient_field_3d * target);
791  N_array_3d * weight_x,
792  N_array_3d * weight_y,
793  N_array_3d * weight_z,
794  N_geom_data * geom,
796  gradfield);
798  field, N_array_3d * x_comp,
799  N_array_3d * y_comp,
800  N_array_3d * z_comp);
801 
802 #endif
void N_free_array_2d(N_array_2d *data)
Release the memory of a N_array_2d structure.
Definition: N_arrays.c:127
N_gradient_neighbours_x * xt
Definition: N_pde.h:622
N_gradient_neighbours_z * zt
Definition: N_pde.h:630
N_gradient_neighbours_3d * N_alloc_gradient_neighbours_3d(void)
Allocate a N_gradient_neighbours_3d structure.
Definition: N_gradient.c:776
N_gradient_neighbours_z * N_create_gradient_neighbours_z(double NWZ, double NZ, double NEZ, double WZ, double CZ, double EZ, double SWZ, double SZ, double SEZ)
Allocate and initialize a N_gradient_neighbours_z structure.
Definition: N_gradient.c:524
void N_free_gradient_neighbours_z(N_gradient_neighbours_z *grad)
Free&#39;s a N_gradient_neighbours_z structure.
Definition: N_gradient.c:500
N_les * N_alloc_nquad_les_Ax(int cols, int rows, int type)
Allocate memory for a (not) quadratic linear equation system which includes the Matrix A and vector x...
Definition: N_les.c:51
void N_free_gradient_2d(N_gradient_2d *grad)
Free&#39;s a N_gradient_2d structure.
Definition: N_gradient.c:42
Matrix entries for a mass balance 5/7/9 star system.
Definition: N_pde.h:342
Gradient between the cell neighbours in X and Y direction.
Definition: N_pde.h:609
float b
Definition: named_colr.c:8
void N_free_gradient_field_3d(N_gradient_field_3d *field)
Free&#39;s a N_gradient_neighbours_3d structure.
Definition: N_gradient.c:1053
int N_is_array_3d_value_null(N_array_3d *data, int col, int row, int depth)
This function returns 1 if value of N_array_3d data at position col, row, depth is of type null...
Definition: N_arrays.c:878
void N_get_array_3d_value(N_array_3d *data, int col, int row, int depth, void *value)
This function writes the value of N_array_3d data at position col, row, depth to the variable value...
Definition: N_arrays.c:826
N_gradient_field_3d * N_compute_gradient_field_3d(N_array_3d *pot, N_array_3d *weight_x, N_array_3d *weight_y, N_array_3d *weight_z, N_geom_data *geom, N_gradient_field_3d *gradfield)
This function computes the gradient based on the input N_array_3d pot (that means potential)...
N_data_star * N_alloc_27star(void)
allocate a 27 point star data structure
int rows
Definition: N_pde.h:162
N_les_callback_3d * N_alloc_les_callback_3d(void)
Allocate the structure holding the callback function.
double N_calc_harmonic_mean_n(double *a, int size)
Calculate the harmonical mean of the values in vector a of size n.
Definition: N_tools.c:142
Gradient between the cell neighbours in Z direction.
Definition: N_pde.h:601
int quad
Definition: N_pde.h:94
double N_calc_geom_mean(double a, double b)
Calculate the geometrical mean of values a and b.
Definition: N_tools.c:77
N_les * N_assemble_les_3d_active(int les_type, N_geom_data *geom, N_array_3d *status, N_array_3d *start_val, void *data, N_les_callback_3d *call)
Assemble a linear equation system (les) based on 3d location data (g3d) active cells.
N_les * N_alloc_nquad_les(int cols, int rows, int type)
Allocate memory for a (not) quadratic linear equation system which includes the Matrix A...
Definition: N_les.c:35
N_array_2d * N_read_rast_to_array_2d(char *name, N_array_2d *array)
Read a raster map into a N_array_2d structure.
Definition: N_arrays_io.c:46
void N_print_gradient_field_3d_info(N_gradient_field_3d *field)
Print gradient field information to stdout.
Definition: N_gradient.c:1098
int dim
Definition: N_pde.h:131
string name
Definition: render.py:1314
void N_free_gradient_neighbours_x(N_gradient_neighbours_x *grad)
Free&#39;s a N_gradient_neighbours_x structure.
Definition: N_gradient.c:309
N_les * N_alloc_les_Ax_b(int rows, int type)
Allocate memory for a quadratic linear equation system which includes the Matrix A, vector x and vector b.
Definition: N_les.c:145
void N_free_array_3d(N_array_3d *data)
Release the memory of a N_array_3d.
Definition: N_arrays.c:777
N_les * N_assemble_les_2d(int les_type, N_geom_data *geom, N_array_2d *status, N_array_2d *start_val, void *data, N_les_callback_2d *call)
Assemble a linear equation system (les) based on 2d location data (raster) and active cells...
#define min(x, y)
Definition: draw2.c:68
void N_compute_gradient_field_components_2d(N_gradient_field_2d *field, N_array_2d *x_comp, N_array_2d *y_comp)
Calculate the x and y vector components from a gradient field for each cell and stores them in the pr...
Gradient between the cells in X and Y direction.
Definition: N_pde.h:525
N_gradient_neighbours_y * y
Definition: N_pde.h:613
void N_put_array_3d_value_null(N_array_3d *data, int col, int row, int depth)
This function writes a null value to the N_array_3d data at position col, row, depth.
Definition: N_arrays.c:1073
N_les * N_alloc_les(int rows, int type)
Allocate memory for a quadratic linear equation system which includes the Matrix A, vector x and vector b.
Definition: N_les.c:100
int N_copy_gradient_3d(N_gradient_3d *source, N_gradient_3d *target)
copy a N_gradient_3d structure
Definition: N_gradient.c:215
double N_calc_geom_mean_n(double *a, int size)
Calculate the geometrical mean of the values in vector a of size n.
Definition: N_tools.c:97
struct Option * N_define_standard_option(int opt)
Create standardised Option structure related to the gpde library.
N_STD_OPT
Standard options of the gpde library.
Definition: N_pde.h:453
FCELL N_get_array_2d_f_value(N_array_2d *data, int col, int row)
Returns the value of type FCELL at position col, row.
Definition: N_arrays.c:343
void N_put_array_2d_d_value(N_array_2d *data, int col, int row, DCELL value)
Writes a DCELL value to the N_array_2d struct at position col, row.
Definition: N_arrays.c:581
void N_put_array_2d_f_value(N_array_2d *data, int col, int row, FCELL value)
Writes a FCELL value to the N_array_2d struct at position col, row.
Definition: N_arrays.c:551
N_gradient_field_2d * N_compute_gradient_field_2d(N_array_2d *pot, N_array_2d *weight_x, N_array_2d *weight_y, N_geom_data *geom, N_gradient_field_2d *gradfield)
This function computes the gradient based on the input N_array_2d pot (potential), a weighting factor N_array_2d named weight and the distance between two cells saved in the N_geom_data struct.
N_array_2d * N_math_array_2d(N_array_2d *a, N_array_2d *b, N_array_2d *result, int type)
Performe calculations with two input arrays, the result is written to a third array.
#define NE
Definition: dataquad.h:18
void N_calc_gradient_field_3d_stats(N_gradient_field_3d *field)
Calculate basic statistics of a gradient field.
int N_copy_gradient_neighbours_2d(N_gradient_neighbours_2d *source, N_gradient_neighbours_2d *target)
copy a N_gradient_neighbours_2d structure
Definition: N_gradient.c:666
void N_calc_array_2d_stats(N_array_2d *a, double *min, double *max, double *sum, int *nonull, int withoffset)
Calculate basic statistics of the N_array_2d struct.
N_array_2d * N_alloc_array_2d(int cols, int rows, int offset, int type)
Allocate memory for a N_array_2d data structure.
Definition: N_arrays.c:69
double WC
Definition: N_pde.h:536
double N_calc_quad_mean(double a, double b)
Calculate the quadratic mean of values a and b.
Definition: N_tools.c:172
#define NW
Definition: dataquad.h:17
double * x
Definition: N_pde.h:88
int y
Definition: plot.c:34
callback structure for 2d matrix assembling
Definition: N_pde.h:364
int depths
Definition: N_pde.h:139
int N_copy_gradient_neighbours_x(N_gradient_neighbours_x *source, N_gradient_neighbours_x *target)
copy a N_gradient_neighbours_x structure
Definition: N_gradient.c:362
N_gradient_3d * N_alloc_gradient_3d(void)
Allocate a N_gradient_3d structure.
Definition: N_gradient.c:152
#define max(x, y)
Definition: draw2.c:69
int N_copy_gradient_2d(N_gradient_2d *source, N_gradient_2d *target)
copy a N_gradient_2d structure
Definition: N_gradient.c:85
N_array_3d * y_array
Definition: N_pde.h:653
N_gradient_2d * N_create_gradient_2d(double NC, double SC, double WC, double EC)
allocate and initialize a N_gradient_2d structure
Definition: N_gradient.c:60
void N_free_gradient_field_2d(N_gradient_field_2d *field)
Free&#39;s a N_gradient_neighbours_2d structure.
Definition: N_gradient.c:952
int offset
Definition: N_pde.h:220
int N_les_pivot_create(N_les *les)
N_array_3d * N_alloc_array_3d(int cols, int rows, int depths, int offset, int type)
Allocate memory for a N_array_3d data structure.
Definition: N_arrays.c:723
N_gradient_neighbours_2d * N_create_gradient_neighbours_2d(N_gradient_neighbours_x *x, N_gradient_neighbours_y *y)
Allocate and initialize a N_gradient_neighbours_2d structure.
Definition: N_gradient.c:632
int count
Definition: N_pde.h:345
double dy
Definition: N_pde.h:134
N_data_star * N_alloc_5star(void)
allocate a 5 point star data structure
N_gradient_neighbours_z * zb
Definition: N_pde.h:631
int offset
Definition: N_pde.h:164
N_gradient_neighbours_2d * N_get_gradient_neighbours_2d(N_gradient_field_2d *field, N_gradient_neighbours_2d *gradient, int col, int row)
Return a N_gradient_neighbours_2d structure calculated from the input gradient field at position [row...
Definition: N_gradient.c:705
tuple size
value.Bind(wx.EVT_TEXT, self.OnVolumeIsosurfMap)
Definition: tools.py:2334
double W_B
Definition: N_pde.h:350
void N_print_array_2d_info(N_array_2d *data)
This function writes the data info of the array data to stdout.
Definition: N_arrays.c:608
CELL * cell_array
Definition: N_pde.h:165
N_data_star * N_alloc_9star(void)
allocate a 9 point star data structure
#define SW
Definition: dataquad.h:19
N_gradient_3d * N_get_gradient_3d(N_gradient_field_3d *field, N_gradient_3d *gradient, int col, int row, int depth)
Return a N_gradient_3d structure calculated from the input gradient field at position [depth][row][co...
Definition: N_gradient.c:249
N_geom_data * N_init_geom_data_2d(struct Cell_head *region, N_geom_data *geodata)
Initiate a pde geometry data structure with a 2d region.
Definition: N_geom.c:115
#define D
Definition: gis/intersect.c:74
#define C
Definition: intr_char.c:17
Geometric information about the structured grid.
Definition: N_pde.h:127
N_les * N_alloc_nquad_les_A(int cols, int rows, int type)
Allocate memory for a (not) quadratic linear equation system which includes the Matrix A...
Definition: N_les.c:67
N_array_2d * y_array
Definition: N_pde.h:641
double dz
Definition: N_pde.h:135
double N_norm_array_3d(N_array_3d *a, N_array_3d *b, int type)
Calculate the norm of the two input arrays.
N_gradient_field_2d * N_alloc_gradient_field_2d(int cols, int rows)
Allocate a N_gradient_field_2d.
Definition: N_gradient.c:928
void N_compute_gradient_field_components_3d(N_gradient_field_3d *field, N_array_3d *x_comp, N_array_3d *y_comp, N_array_3d *z_comp)
Calculate the x, y and z vector components from a gradient field for each cell and store them in the ...
tuple data
int rows
Definition: N_pde.h:140
void N_set_les_callback_3d_func(N_les_callback_3d *data, N_data_star *(*callback_func_3d)())
Set the callback function which is called while assembling the les in 3d.
N_gradient_neighbours_x * xb
Definition: N_pde.h:624
N_data_star * N_create_7star(double C, double W, double E, double N, double S, double T, double B, double V)
allocate and initialize a 7 point star data structure
double ** A
Definition: N_pde.h:90
N_data_star * N_create_27star(double C, double W, double E, double N, double S, double NW, double SW, double NE, double SE, double T, double W_T, double E_T, double N_T, double S_T, double NW_T, double SW_T, double NE_T, double SE_T, double B, double W_B, double E_B, double N_B, double S_B, double NW_B, double SW_B, double NE_B, double SE_B, double V)
allocate and initialize a 27 point star data structure
int N_convert_array_3d_null_to_zero(N_array_3d *a)
Convert all null values to zero values.
N_gradient_neighbours_y * yb
Definition: N_pde.h:628
void N_print_array_3d_info(N_array_3d *data)
Write the info of the array to stdout.
Definition: N_arrays.c:1194
Gradient between the cell neighbours in X, Y and Z direction.
Definition: N_pde.h:619
N_array_3d * z_array
Definition: N_pde.h:654
Gradient between the cell neighbours in X direction.
Definition: N_pde.h:585
double W_T
Definition: N_pde.h:348
double N_calc_quad_mean_n(double *a, int size)
Calculate the quadratic mean of the values in vector a of size n.
Definition: N_tools.c:192
double N_exp_upwinding(double sprod, double distance, double D)
exponential upwinding stabilization algorithm
Definition: N_upwind.c:63
N_gradient_neighbours_y * yt
Definition: N_pde.h:626
FCELL * fcell_array
Definition: N_pde.h:166
N_gradient_neighbours_x * N_create_gradient_neighbours_x(double NWN, double NEN, double WC, double EC, double SWS, double SES)
Allocate and initialize a N_gradient_neighbours_x structure.
Definition: N_gradient.c:331
N_gradient_neighbours_z * N_alloc_gradient_neighbours_z(void)
Allocate a N_gradient_neighbours_z structure.
Definition: N_gradient.c:483
N_les * N_assemble_les_2d_param(int les_type, N_geom_data *geom, N_array_2d *status, N_array_2d *start_val, void *data, N_les_callback_2d *call, int cell_type)
Assemble a linear equation system (les) based on 2d location data (raster)
double Az
Definition: N_pde.h:137
N_gradient_3d * N_create_gradient_3d(double NC, double SC, double WC, double EC, double TC, double BC)
allocate and initialize a N_gradient_3d structure
Definition: N_gradient.c:188
int N_convert_array_2d_null_to_zero(N_array_2d *a)
Convert all null values to zero values.
G_math_spvector ** Asp
Definition: N_pde.h:91
void N_put_array_2d_value_null(N_array_2d *data, int col, int row)
Writes the null value to the N_array_2d struct at position col, row.
Definition: N_arrays.c:455
void N_print_gradient_field_2d_info(N_gradient_field_2d *field)
Print gradient field information to stdout.
Definition: N_gradient.c:994
void N_copy_array_2d(N_array_2d *source, N_array_2d *target)
Copy the source N_array_2d struct to the target N_array_2d struct.
Definition: N_arrays_calc.c:43
void N_put_array_2d_c_value(N_array_2d *data, int col, int row, CELL value)
Writes a CELL value to the N_array_2d struct at position col, row.
Definition: N_arrays.c:521
N_data_star * N_create_9star(double C, double W, double E, double N, double S, double NW, double SW, double NE, double SE, double V)
allocate and initialize a 9 point star data structure
N_gradient_neighbours_2d * N_alloc_gradient_neighbours_2d(void)
Allocate a N_gradient_neighbours_2d structure.
Definition: N_gradient.c:589
N_les * N_assemble_les_2d_dirichlet(int les_type, N_geom_data *geom, N_array_2d *status, N_array_2d *start_val, void *data, N_les_callback_2d *call)
Assemble a linear equation system (les) based on 2d location data (raster) and active and dirichlet c...
char * value
Definition: env.c:30
float N_get_array_3d_f_value(N_array_3d *data, int col, int row, int depth)
This function returns the value of type float at position col, row, depth.
Definition: N_arrays.c:958
DCELL * dcell_array
Definition: N_pde.h:167
int N_les_integrate_dirichlet_3d(N_les *les, N_geom_data *geom, N_array_3d *status, N_array_3d *start_val)
Integrate Dirichlet or Transmission boundary conditions into the les (3d)
N_array_3d * N_read_rast3d_to_array_3d(char *name, N_array_3d *array, int mask)
Read a volume map into a N_array_3d structure.
Definition: N_arrays_io.c:276
N_les * N_assemble_les_3d(int les_type, N_geom_data *geom, N_array_3d *status, N_array_3d *start_val, void *data, N_les_callback_3d *call)
Assemble a linear equation system (les) based on 3d location data (g3d) active cells.
void N_print_array_2d(N_array_2d *data)
Write info and content of the N_array_2d struct to stdout.
Definition: N_arrays.c:634
N_data_star * N_create_5star(double C, double W, double E, double N, double S, double V)
allocate and initialize a 5 point star data structure
N_gradient_neighbours_3d * N_create_gradient_neighbours_3d(N_gradient_neighbours_x *xt, N_gradient_neighbours_x *xc, N_gradient_neighbours_x *xb, N_gradient_neighbours_y *yt, N_gradient_neighbours_y *yc, N_gradient_neighbours_y *yb, N_gradient_neighbours_z *zt, N_gradient_neighbours_z *zb)
Allocate and initialize a N_gradient_neighbours_3d structure.
Definition: N_gradient.c:832
int cols
Definition: N_pde.h:141
Gradient between the cells in X, Y and Z direction.
Definition: N_pde.h:533
N_data_star * N_callback_template_3d(void *data, N_geom_data *geom, int col, int row, int depth)
A callback template creates a 7 point star structure.
void N_put_array_3d_f_value(N_array_3d *data, int col, int row, int depth, float value)
This function writes a float value to the N_array_3d data at position col, row, depth.
Definition: N_arrays.c:1145
void N_print_array_3d(N_array_3d *data)
Write info and content of the array data to stdout.
Definition: N_arrays.c:1220
void N_calc_array_3d_stats(N_array_3d *a, double *min, double *max, double *sum, int *nonull, int withoffset)
Calculate basic statistics of the N_array_3d struct.
double N_get_geom_data_area_of_cell(N_geom_data *geom, int row)
Get the areay size in square meter of one cell (x*y) at row.
Definition: N_geom.c:199
int type
Definition: N_pde.h:344
#define SE
Definition: dataquad.h:20
double N_calc_arith_mean(double a, double b)
Calculate the arithmetic mean of values a and b.
Definition: N_tools.c:33
float * fcell_array
Definition: N_pde.h:221
void N_get_array_2d_value(N_array_2d *data, int col, int row, void *value)
Write the value of the N_array_2d struct at position col, row to value.
Definition: N_arrays.c:178
void N_set_les_callback_2d_func(N_les_callback_2d *data, N_data_star *(*callback_func_2d)())
Set the callback function which is called while assembling the les in 2d.
N_les * N_assemble_les_3d_param(int les_type, N_geom_data *geom, N_array_3d *status, N_array_3d *start_val, void *data, N_les_callback_3d *call, int cell_type)
Assemble a linear equation system (les) based on 3d location data (g3d)
N_data_star * N_alloc_7star(void)
allocate a 7 point star data structure
N_array_3d * N_math_array_3d(N_array_3d *a, N_array_3d *b, N_array_3d *result, int type)
Performe calculations with two input arrays, the result is written to a third array.
N_les * N_alloc_les_param(int cols, int rows, int type, int parts)
Allocate memory for a quadratic or not quadratic linear equation system.
Definition: N_les.c:178
int planimetric
Definition: N_pde.h:129
void N_write_array_3d_to_rast3d(N_array_3d *array, char *name, int mask)
Write a N_array_3d struct to a volume map.
Definition: N_arrays_io.c:408
int N_is_array_2d_value_null(N_array_2d *data, int col, int row)
Returns 1 if the value of N_array_2d struct at postion col, row is of type null, otherwise 0...
Definition: N_arrays.c:228
void N_free_gradient_neighbours_2d(N_gradient_neighbours_2d *grad)
Free&#39;s a N_gradient_neighbours_2d structure.
Definition: N_gradient.c:610
N_geom_data * N_alloc_geom_data(void)
Allocate the pde geometry data structure and return a pointer to the new allocated structure...
Definition: N_geom.c:30
double * dcell_array
Definition: N_pde.h:222
N_gradient_neighbours_x * N_alloc_gradient_neighbours_x(void)
Allocate a N_gradient_neighbours_x structure.
Definition: N_gradient.c:292
double N_get_array_3d_d_value(N_array_3d *data, int col, int row, int depth)
This function returns the value of type float at position col, row, depth.
Definition: N_arrays.c:987
N_gradient_2d * N_get_gradient_2d(N_gradient_field_2d *field, N_gradient_2d *gradient, int col, int row)
Return a N_gradient_2d structure calculated from the input gradient field at position [row][col]...
Definition: N_gradient.c:115
int type
Definition: N_pde.h:161
N_les * N_alloc_les_A(int rows, int type)
Allocate memory for a quadratic linear equation system which includes the Matrix A.
Definition: N_les.c:130
int cols
Definition: N_pde.h:93
tuple cols
void N_free_geom_data(N_geom_data *geom)
Release memory of a pde geometry data structure.
Definition: N_geom.c:50
void N_put_array_2d_value(N_array_2d *data, int col, int row, char *value)
Writes a value to the N_array_2d struct at position col, row.
Definition: N_arrays.c:408
N_les * N_alloc_nquad_les_Ax_b(int cols, int rows, int type)
Allocate memory for a (not) quadratic linear equation system which includes the Matrix A...
Definition: N_les.c:83
int N_copy_gradient_neighbours_z(N_gradient_neighbours_z *source, N_gradient_neighbours_z *target)
copy a N_gradient_neighbours_z structure
Definition: N_gradient.c:559
double * area
Definition: N_pde.h:130
N_array_3d * x_array
Definition: N_pde.h:652
N_les * N_alloc_les_Ax(int rows, int type)
Allocate memory for a quadratic linear equation system which includes the Matrix A and vector x...
Definition: N_les.c:115
N_gradient_field_3d * N_alloc_gradient_field_3d(int cols, int rows, int depths)
Allocate a N_gradient_field_3d.
Definition: N_gradient.c:1026
double WC
Definition: N_pde.h:528
DCELL N_get_array_2d_d_value(N_array_2d *data, int col, int row)
Returns the value of type DCELL at position col, row.
Definition: N_arrays.c:375
void N_write_array_2d_to_rast(N_array_2d *array, char *name)
Write a N_array_2d struct to a raster map.
Definition: N_arrays_io.c:181
N_data_star * N_callback_template_2d(void *data, N_geom_data *geom, int col, int row)
A callback template creates a 9 point star structure.
void N_free_gradient_3d(N_gradient_3d *grad)
Free&#39;s a N_gradient_3d structure.
Definition: N_gradient.c:167
void N_put_array_3d_d_value(N_array_3d *data, int col, int row, int depth, double value)
Writes a double value to the N_array_3d struct at position col, row, depth.
Definition: N_arrays.c:1172
N_gradient_neighbours_x * xc
Definition: N_pde.h:623
int N_get_array_3d_type(N_array_3d *array)
Return the data type of the N_array_3d.
Definition: N_arrays.c:806
int N_get_array_2d_type(N_array_2d *array)
Return the data type of the N_array_2d struct.
Definition: N_arrays.c:161
#define N
Definition: inverse.c:8
N_gradient_neighbours_y * N_create_gradient_neighbours_y(double NWW, double NEE, double NC, double SC, double SWW, double SEE)
Allocate and initialize a N_gradient_neighbours_y structure.
Definition: N_gradient.c:426
double * b
Definition: N_pde.h:89
double N_norm_array_2d(N_array_2d *a, N_array_2d *b, int type)
Calculate the norm of the two input arrays.
CELL N_get_array_2d_c_value(N_array_2d *data, int col, int row)
Returns the value of type CELL at position col, row.
Definition: N_arrays.c:311
N_gradient_neighbours_y * yc
Definition: N_pde.h:627
double N_calc_harmonic_mean(double a, double b)
Calculate the harmonical mean of values a and b.
Definition: N_tools.c:121
int N_copy_gradient_neighbours_3d(N_gradient_neighbours_3d *source, N_gradient_neighbours_3d *target)
copy a N_gradient_neighbours_3d structure
Definition: N_gradient.c:883
void N_free_gradient_neighbours_3d(N_gradient_neighbours_3d *grad)
Free&#39;s a N_gradient_neighbours_3d structure.
Definition: N_gradient.c:803
N_gradient_2d * N_alloc_gradient_2d(void)
Allocate a N_gradient_2d structure.
Definition: N_gradient.c:27
callback structure for 3d matrix assembling
Definition: N_pde.h:356
Gradient between the cell neighbours in Y direction.
Definition: N_pde.h:593
void N_calc_gradient_field_2d_stats(N_gradient_field_2d *field)
Calculate basic statistics of a gradient field.
double W
Definition: N_pde.h:346
int N_les_integrate_dirichlet_2d(N_les *les, N_geom_data *geom, N_array_2d *status, N_array_2d *start_val)
Integrate Dirichlet or Transmission boundary conditions into the les (2s)
double dx
Definition: N_pde.h:133
int rows
Definition: N_pde.h:92
N_gradient_neighbours_x * x
Definition: N_pde.h:612
N_les * N_assemble_les_3d_dirichlet(int les_type, N_geom_data *geom, N_array_3d *status, N_array_3d *start_val, void *data, N_les_callback_3d *call)
Assemble a linear equation system (les) based on 3d location data (g3d) active and dirichlet cells...
void N_put_array_3d_value(N_array_3d *data, int col, int row, int depth, char *value)
This function writes a value to the N_array_3d data at position col, row, depth.
Definition: N_arrays.c:1018
int N_copy_gradient_field_2d(N_gradient_field_2d *source, N_gradient_field_2d *target)
Copy N_gradient_field_2d structure from source to target.
Definition: N_gradient.c:974
int type
Definition: N_pde.h:217
int rows_intern
Definition: N_pde.h:219
int N_copy_gradient_field_3d(N_gradient_field_3d *source, N_gradient_field_3d *target)
Copy N_gradient_field_3d structure from source to target.
Definition: N_gradient.c:1077
SC
Definition: wxgui.py:32
N_les * N_assemble_les_2d_active(int les_type, N_geom_data *geom, N_array_2d *status, N_array_2d *start_val, void *data, N_les_callback_2d *call)
Assemble a linear equation system (les) based on 2d location data (raster) and active cells...
void N_free_les(N_les *les)
Release the memory of the linear equation system.
Definition: N_les.c:304
The linear equation system (les) structure.
Definition: N_pde.h:86
N_geom_data * N_init_geom_data_3d(G3D_Region *region3d, N_geom_data *geodata)
Initiate a pde geometry data structure with a 3d region.
Definition: N_geom.c:73
N_gradient_neighbours_y * N_alloc_gradient_neighbours_y(void)
Allocate a N_gradient_neighbours_y structure.
Definition: N_gradient.c:388
N_array_2d * x_array
Definition: N_pde.h:640
int rows
Definition: N_pde.h:218
int type
Definition: N_pde.h:95
void N_copy_array_3d(N_array_3d *source, N_array_3d *target)
Copy the source N_array_3d struct to the target N_array_3d struct.
double N_full_upwinding(double sprod, double distance, double D)
full upwinding stabilization algorithm
Definition: N_upwind.c:33
int rows_intern
Definition: N_pde.h:163
void N_print_les(N_les *les)
prints the linear equation system to stdout
Definition: N_les.c:252
int N_copy_gradient_neighbours_y(N_gradient_neighbours_y *source, N_gradient_neighbours_y *target)
copy a N_gradient_neighbours_y structure
Definition: N_gradient.c:457
N_les_callback_2d * N_alloc_les_callback_2d(void)
Allocate the structure holding the callback function.
void N_free_gradient_neighbours_y(N_gradient_neighbours_y *grad)
Free&#39;s a N_gradient_neighbours_y structure.
Definition: N_gradient.c:405
double N_calc_arith_mean_n(double *a, int size)
Calculate the arithmetic mean of the values in vector a of size n.
Definition: N_tools.c:53