GRASS Programmer's Manual  6.5.svn(2012)-r51648
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines
N_solute_transport.c
Go to the documentation of this file.
00001 
00002 /*****************************************************************************
00003 *
00004 * MODULE:       Grass PDE Numerical Library
00005 * AUTHOR(S):    Soeren Gebbert, Berlin (GER) Dec 2006
00006 *               soerengebbert <at> gmx <dot> de
00007 *               
00008 * PURPOSE:      solute transport in porous media
00009 *               part of the gpde library
00010 *
00011 * COPYRIGHT:    (C) 2007 by the GRASS Development Team
00012 *
00013 *               This program is free software under the GNU General Public
00014 *               License (>=v2). Read the file COPYING that comes with GRASS
00015 *               for details.
00016 *
00017 *****************************************************************************/
00018 
00019 #include "grass/N_solute_transport.h"
00020 #include <math.h>
00021 
00022 /* ************************************************************************* *
00023  * ************************************************************************* *
00024  * ************************************************************************* */
00028 N_data_star *N_callback_solute_transport_3d(void *solutedata,
00029                                             N_geom_data * geom, int col,
00030                                             int row, int depth)
00031 {
00032     double Df_e = 0, Df_w = 0, Df_n = 0, Df_s = 0, Df_t = 0, Df_b = 0;
00033 
00034     double dx, dy, dz, Az;
00035 
00036     double diff_x, diff_y, diff_z;
00037 
00038     double diff_xw, diff_yn;
00039 
00040     double diff_xe, diff_ys;
00041 
00042     double diff_zt, diff_zb;
00043 
00044     double cin = 0, cg, cg_start;
00045 
00046     double R, nf, cs, q;
00047 
00048     double C, W, E, N, S, T, B, V;
00049 
00050     double vw = 0, ve = 0, vn = 0, vs = 0, vt = 0, vb = 0;
00051 
00052     double Ds_w = 0, Ds_e = 0, Ds_n = 0, Ds_s = 0, Ds_t = 0, Ds_b = 0;
00053 
00054     double Dw = 0, De = 0, Dn = 0, Ds = 0, Dt = 0, Db = 0;
00055 
00056     double rw = 0.5, re = 0.5, rn = 0.5, rs = 0.5, rt = 0.5, rb = 0.5;
00057 
00058     N_solute_transport_data3d *data = NULL;
00059 
00060     N_data_star *mat_pos;
00061 
00062     N_gradient_3d grad;
00063 
00064     /*cast the void pointer to the right data structure */
00065     data = (N_solute_transport_data3d *) solutedata;
00066 
00067     N_get_gradient_3d(data->grad, &grad, col, row, depth);
00068 
00069     dx = geom->dx;
00070     dy = geom->dy;
00071     dz = geom->dz;
00072     Az = N_get_geom_data_area_of_cell(geom, row);
00073 
00074     /*read the data from the arrays */
00075     cg_start = N_get_array_3d_d_value(data->c_start, col, row, depth);
00076     cg = N_get_array_3d_d_value(data->c, col, row, depth);
00077 
00078     /*get the surrounding diffusion tensor entries */
00079     diff_x = N_get_array_3d_d_value(data->diff_x, col, row, depth);
00080     diff_y = N_get_array_3d_d_value(data->diff_y, col, row, depth);
00081     diff_z = N_get_array_3d_d_value(data->diff_z, col, row, depth);
00082     diff_xw = N_get_array_3d_d_value(data->diff_x, col - 1, row, depth);
00083     diff_xe = N_get_array_3d_d_value(data->diff_x, col + 1, row, depth);
00084     diff_yn = N_get_array_3d_d_value(data->diff_y, col, row - 1, depth);
00085     diff_ys = N_get_array_3d_d_value(data->diff_y, col, row + 1, depth);
00086     diff_zt = N_get_array_3d_d_value(data->diff_z, col, row, depth + 1);
00087     diff_zb = N_get_array_3d_d_value(data->diff_z, col, row, depth - 1);
00088 
00089     /* calculate the diffusion on the cell borders using the harmonical mean */
00090     Df_w = N_calc_harmonic_mean(diff_xw, diff_x);
00091     Df_e = N_calc_harmonic_mean(diff_xe, diff_x);
00092     Df_n = N_calc_harmonic_mean(diff_yn, diff_y);
00093     Df_s = N_calc_harmonic_mean(diff_ys, diff_y);
00094     Df_t = N_calc_harmonic_mean(diff_zt, diff_z);
00095     Df_b = N_calc_harmonic_mean(diff_zb, diff_z);
00096 
00097     /* calculate the dispersion */
00098     /*todo */
00099 
00100     /* calculate the velocity parts  with full upwinding scheme */
00101     vw = grad.WC;
00102     ve = grad.EC;
00103     vn = grad.NC;
00104     vs = grad.SC;
00105     vt = grad.TC;
00106     vb = grad.BC;
00107 
00108     /* put the diffusion and dispersion together */
00109     Dw = ((Df_w + Ds_w)) / dx;
00110     De = ((Df_e + Ds_e)) / dx;
00111     Dn = ((Df_n + Ds_n)) / dy;
00112     Ds = ((Df_s + Ds_s)) / dy;
00113     Dt = ((Df_t + Ds_t)) / dz;
00114     Db = ((Df_b + Ds_b)) / dz;
00115 
00116     rw = N_exp_upwinding(-1 * vw, dx, Dw);
00117     re = N_exp_upwinding(ve, dx, De);
00118     rs = N_exp_upwinding(-1 * vs, dy, Ds);
00119     rn = N_exp_upwinding(vn, dy, Dn);
00120     rb = N_exp_upwinding(-1 * vb, dz, Dn);
00121     rt = N_exp_upwinding(vt, dz, Dn);
00122 
00123     /*mass balance center cell to western cell */
00124     W = -1 * (Dw) * dy * dz - vw * (1 - rw) * dy * dz;
00125     /*mass balance center cell to eastern cell */
00126     E = -1 * (De) * dy * dz + ve * (1 - re) * dy * dz;
00127     /*mass balance center cell to southern cell */
00128     S = -1 * (Ds) * dx * dz - vs * (1 - rs) * dx * dz;
00129     /*mass balance center cell to northern cell */
00130     N = -1 * (Dn) * dx * dz + vn * (1 - rn) * dx * dz;
00131     /*mass balance center cell to bottom cell */
00132     B = -1 * (Db) * Az - vb * (1 - rb) * Az;
00133     /*mass balance center cell to top cell */
00134     T = -1 * (Dt) * Az + vt * (1 - rt) * Az;
00135 
00136     /* Retardation */
00137     R = N_get_array_3d_d_value(data->R, col, row, depth);
00138     /* Inner sources */
00139     cs = N_get_array_3d_d_value(data->cs, col, row, depth);
00140     /* effective porosity */
00141     nf = N_get_array_3d_d_value(data->nf, col, row, depth);
00142     /* groundwater sources and sinks */
00143     q = N_get_array_3d_d_value(data->q, col, row, depth);
00144     /* concentration of influent water */
00145     cin = N_get_array_3d_d_value(data->cin, col, row, depth);
00146 
00147     /*the diagonal entry of the matrix */
00148     C = ((Dw - vw) * dy * dz +
00149          (De + ve) * dy * dz +
00150          (Ds - vs) * dx * dz +
00151          (Dn + vn) * dx * dz +
00152          (Db - vb) * Az + (Dt + vt) * Az + Az * dz * R / data->dt - q / nf);
00153 
00154     /*the entry in the right side b of Ax = b */
00155     V = (cs + cg_start * Az * dz * R / data->dt - q / nf * cin);
00156 
00157     /*
00158      * printf("nf %g\n", nf);
00159      * printf("q %g\n", q);
00160      * printf("cs %g\n", cs);
00161      * printf("cin %g\n", cin);
00162      * printf("cg %g\n", cg);
00163      * printf("cg_start %g\n", cg_start);
00164      * printf("Az %g\n", Az);
00165      * printf("z %g\n", z);
00166      * printf("R %g\n", R);
00167      * printf("dt %g\n", data->dt);
00168      */
00169     G_debug(6, "N_callback_solute_transport_3d: called [%i][%i][%i]", row,
00170             col, depth);
00171 
00172     /*create the 7 point star entries */
00173     mat_pos = N_create_7star(C, W, E, N, S, T, B, V);
00174 
00175     return mat_pos;
00176 }
00177 
00178 /* ************************************************************************* *
00179  * ************************************************************************* *
00180  * ************************************************************************* */
00198 N_data_star *N_callback_solute_transport_2d(void *solutedata,
00199                                             N_geom_data * geom, int col,
00200                                             int row)
00201 {
00202     double Df_e = 0, Df_w = 0, Df_n = 0, Df_s = 0;
00203 
00204     double z_e = 0, z_w = 0, z_n = 0, z_s = 0;
00205 
00206     double dx, dy, Az;
00207 
00208     double diff_x, diff_y;
00209 
00210     double disp_x, disp_y;
00211 
00212     double z;
00213 
00214     double diff_xw, diff_yn;
00215 
00216     double disp_xw, disp_yn;
00217 
00218     double z_xw, z_yn;
00219 
00220     double diff_xe, diff_ys;
00221 
00222     double disp_xe, disp_ys;
00223 
00224     double z_xe, z_ys;
00225 
00226     double cin = 0, cg, cg_start;
00227 
00228     double R, nf, cs, q;
00229 
00230     double C, W, E, N, S, V, NE, NW, SW, SE;
00231 
00232     double vw = 0, ve = 0, vn = 0, vs = 0;
00233 
00234     double Ds_w = 0, Ds_e = 0, Ds_n = 0, Ds_s = 0;
00235 
00236     double Dw = 0, De = 0, Dn = 0, Ds = 0;
00237 
00238     double rw = 0.5, re = 0.5, rn = 0.5, rs = 0.5;
00239 
00240     N_solute_transport_data2d *data = NULL;
00241 
00242     N_data_star *mat_pos;
00243 
00244     N_gradient_2d grad;
00245 
00246     /*cast the void pointer to the right data structure */
00247     data = (N_solute_transport_data2d *) solutedata;
00248 
00249     N_get_gradient_2d(data->grad, &grad, col, row);
00250 
00251     dx = geom->dx;
00252     dy = geom->dy;
00253     Az = N_get_geom_data_area_of_cell(geom, row);
00254 
00255     /*read the data from the arrays */
00256     cg_start = N_get_array_2d_d_value(data->c_start, col, row);
00257     cg = N_get_array_2d_d_value(data->c, col, row);
00258 
00259     /* calculate the cell height */
00260     z = N_get_array_2d_d_value(data->top, col,
00261                                row) -
00262         N_get_array_2d_d_value(data->bottom, col, row);
00263     z_xw =
00264         N_get_array_2d_d_value(data->top, col - 1,
00265                                row) -
00266         N_get_array_2d_d_value(data->bottom, col - 1, row);
00267     z_xe =
00268         N_get_array_2d_d_value(data->top, col + 1,
00269                                row) -
00270         N_get_array_2d_d_value(data->bottom, col + 1, row);
00271     z_yn =
00272         N_get_array_2d_d_value(data->top, col,
00273                                row - 1) -
00274         N_get_array_2d_d_value(data->bottom, col, row - 1);
00275     z_ys =
00276         N_get_array_2d_d_value(data->top, col,
00277                                row + 1) -
00278         N_get_array_2d_d_value(data->bottom, col, row + 1);
00279 
00280     /*geometrical mean of cell height */
00281     z_w = N_calc_geom_mean(z_xw, z);
00282     z_e = N_calc_geom_mean(z_xe, z);
00283     z_n = N_calc_geom_mean(z_yn, z);
00284     z_s = N_calc_geom_mean(z_ys, z);
00285 
00286     /*get the surrounding diffusion tensor entries */
00287     diff_x = N_get_array_2d_d_value(data->diff_x, col, row);
00288     diff_y = N_get_array_2d_d_value(data->diff_y, col, row);
00289     diff_xw = N_get_array_2d_d_value(data->diff_x, col - 1, row);
00290     diff_xe = N_get_array_2d_d_value(data->diff_x, col + 1, row);
00291     diff_yn = N_get_array_2d_d_value(data->diff_y, col, row - 1);
00292     diff_ys = N_get_array_2d_d_value(data->diff_y, col, row + 1);
00293 
00294     /* calculate the diffusion at the cell borders using the harmonical mean */
00295     Df_w = N_calc_harmonic_mean(diff_xw, diff_x);
00296     Df_e = N_calc_harmonic_mean(diff_xe, diff_x);
00297     Df_n = N_calc_harmonic_mean(diff_yn, diff_y);
00298     Df_s = N_calc_harmonic_mean(diff_ys, diff_y);
00299 
00300     /* calculate the dispersion */
00301     /*get the surrounding dispersion tensor entries */
00302     disp_x = N_get_array_2d_d_value(data->disp_xx, col, row);
00303     disp_y = N_get_array_2d_d_value(data->disp_yy, col, row);
00304     if (N_get_array_2d_d_value(data->status, col - 1, row) ==
00305         N_CELL_TRANSMISSION) {
00306         disp_xw = disp_x;
00307     }
00308     else {
00309         disp_xw = N_get_array_2d_d_value(data->disp_xx, col - 1, row);
00310     }
00311     if (N_get_array_2d_d_value(data->status, col + 1, row) ==
00312         N_CELL_TRANSMISSION) {
00313         disp_xe = disp_x;
00314     }
00315     else {
00316         disp_xe = N_get_array_2d_d_value(data->disp_xx, col + 1, row);
00317     }
00318     if (N_get_array_2d_d_value(data->status, col, row - 1) ==
00319         N_CELL_TRANSMISSION) {
00320         disp_yn = disp_y;
00321     }
00322     else {
00323         disp_yn = N_get_array_2d_d_value(data->disp_yy, col, row - 1);
00324     }
00325     if (N_get_array_2d_d_value(data->status, col, row + 1) ==
00326         N_CELL_TRANSMISSION) {
00327         disp_ys = disp_y;
00328     }
00329     else {
00330         disp_ys = N_get_array_2d_d_value(data->disp_yy, col, row + 1);
00331     }
00332 
00333     /* calculate the dispersion at the cell borders using the harmonical mean */
00334     Ds_w = N_calc_harmonic_mean(disp_xw, disp_x);
00335     Ds_e = N_calc_harmonic_mean(disp_xe, disp_x);
00336     Ds_n = N_calc_harmonic_mean(disp_yn, disp_y);
00337     Ds_s = N_calc_harmonic_mean(disp_ys, disp_y);
00338 
00339     /* put the diffusion and dispersion together */
00340     Dw = ((Df_w + Ds_w)) / dx;
00341     De = ((Df_e + Ds_e)) / dx;
00342     Ds = ((Df_s + Ds_s)) / dy;
00343     Dn = ((Df_n + Ds_n)) / dy;
00344 
00345     vw = -1.0 * grad.WC;
00346     ve = grad.EC;
00347     vs = -1.0 * grad.SC;
00348     vn = grad.NC;
00349 
00350     if (data->stab == N_UPWIND_FULL) {
00351         rw = N_full_upwinding(vw, dx, Dw);
00352         re = N_full_upwinding(ve, dx, De);
00353         rs = N_full_upwinding(vs, dy, Ds);
00354         rn = N_full_upwinding(vn, dy, Dn);
00355     }
00356     else if (data->stab == N_UPWIND_EXP) {
00357         rw = N_exp_upwinding(vw, dx, Dw);
00358         re = N_exp_upwinding(ve, dx, De);
00359         rs = N_exp_upwinding(vs, dy, Ds);
00360         rn = N_exp_upwinding(vn, dy, Dn);
00361     }
00362 
00363     /*mass balance center cell to western cell */
00364     W = -1 * (Dw) * dy * z_w + vw * (1 - rw) * dy * z_w;
00365     /*mass balance center cell to eastern cell */
00366     E = -1 * (De) * dy * z_e + ve * (1 - re) * dy * z_e;
00367     /*mass balance center cell to southern cell */
00368     S = -1 * (Ds) * dx * z_s + vs * (1 - rs) * dx * z_s;
00369     /*mass balance center cell to northern cell */
00370     N = -1 * (Dn) * dx * z_n + vn * (1 - rn) * dx * z_n;
00371 
00372     NW = 0.0;
00373     SW = 0.0;
00374     NE = 0.0;
00375     SE = 0.0;
00376 
00377     /* Retardation */
00378     R = N_get_array_2d_d_value(data->R, col, row);
00379     /* Inner sources */
00380     cs = N_get_array_2d_d_value(data->cs, col, row);
00381     /* effective porosity */
00382     nf = N_get_array_2d_d_value(data->nf, col, row);
00383     /* groundwater sources and sinks */
00384     q = N_get_array_2d_d_value(data->q, col, row);
00385     /* concentration of influent water */
00386     cin = N_get_array_2d_d_value(data->cin, col, row);
00387 
00388     /*the diagonal entry of the matrix */
00389     C = (Dw + vw * rw) * dy * z_w +
00390         (De + ve * re) * dy * z_e +
00391         (Ds + vs * rs) * dx * z_s +
00392         (Dn + vn * rn) * dx * z_n + Az * z * R / data->dt - q / nf;
00393 
00394     /*the entry in the right side b of Ax = b */
00395     V = (cs + cg_start * Az * z * R / data->dt + q / nf * cin);
00396 
00397     /*
00398        fprintf(stderr, "nf %g\n", nf);
00399        fprintf(stderr, "q %g\n", q);
00400        fprintf(stderr, "cs %g\n", cs);
00401        fprintf(stderr, "cin %g\n", cin);
00402        fprintf(stderr, "cg %g\n", cg);
00403        fprintf(stderr, "cg_start %g\n", cg_start);
00404        fprintf(stderr, "Az %g\n", Az);
00405        fprintf(stderr, "z %g\n", z);
00406        fprintf(stderr, "R %g\n", R);
00407        fprintf(stderr, "dt %g\n", data->dt);
00408      */
00409 
00410     G_debug(6, "N_callback_solute_transport_2d: called [%i][%i]", row, col);
00411 
00412     /*create the 9 point star entries */
00413     mat_pos = N_create_9star(C, W, E, N, S, NW, SW, NE, SE, V);
00414 
00415     return mat_pos;
00416 }
00417 
00418 /* ************************************************************************* *
00419  * ************************************************************************* *
00420  * ************************************************************************* */
00436 N_solute_transport_data3d *N_alloc_solute_transport_data3d(int cols, int rows,
00437                                                            int depths)
00438 {
00439     N_solute_transport_data3d *data = NULL;
00440 
00441     data =
00442         (N_solute_transport_data3d *) G_calloc(1,
00443                                                sizeof
00444                                                (N_solute_transport_data3d));
00445 
00446     data->c = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00447     data->c_start = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00448     data->status = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00449     data->diff_x = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00450     data->diff_y = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00451     data->diff_z = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00452     data->q = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00453     data->cs = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00454     data->R = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00455     data->nf = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00456     data->cin = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00457 
00458     /*Allocate the dispersivity tensor */
00459     data->disp_xx = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00460     data->disp_yy = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00461     data->disp_zz = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00462     data->disp_xy = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00463     data->disp_xz = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00464     data->disp_yz = N_alloc_array_3d(cols, rows, depths, 1, DCELL_TYPE);
00465 
00466 
00467     data->grad = N_alloc_gradient_field_3d(cols, rows, depths);
00468     data->stab = N_UPWIND_EXP;
00469 
00470     return data;
00471 }
00472 
00473 /* ************************************************************************* *
00474  * ************************************************************************* *
00475  * ************************************************************************* */
00491 N_solute_transport_data2d *N_alloc_solute_transport_data2d(int cols, int rows)
00492 {
00493     N_solute_transport_data2d *data = NULL;
00494 
00495     data =
00496         (N_solute_transport_data2d *) G_calloc(1,
00497                                                sizeof
00498                                                (N_solute_transport_data2d));
00499 
00500     data->c = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00501     data->c_start = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00502     data->status = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00503     data->diff_x = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00504     data->diff_y = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00505     data->q = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00506     data->cs = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00507     data->R = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00508     data->nf = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00509     data->cin = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00510     data->top = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00511     data->bottom = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00512 
00513     /*Allocate the dispersivity tensor */
00514     data->disp_xx = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00515     data->disp_yy = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00516     data->disp_xy = N_alloc_array_2d(cols, rows, 1, DCELL_TYPE);
00517 
00518     data->grad = N_alloc_gradient_field_2d(cols, rows);
00519     data->stab = N_UPWIND_EXP;
00520 
00521     return data;
00522 }
00523 
00524 /* ************************************************************************* *
00525  * ************************************************************************* *
00526  * ************************************************************************* */
00533 void N_free_solute_transport_data3d(N_solute_transport_data3d * data)
00534 {
00535     N_free_array_3d(data->c);
00536     N_free_array_3d(data->c_start);
00537     N_free_array_3d(data->status);
00538     N_free_array_3d(data->diff_x);
00539     N_free_array_3d(data->diff_y);
00540     N_free_array_3d(data->diff_z);
00541     N_free_array_3d(data->q);
00542     N_free_array_3d(data->cs);
00543     N_free_array_3d(data->R);
00544     N_free_array_3d(data->nf);
00545     N_free_array_3d(data->cin);
00546 
00547     N_free_array_3d(data->disp_xx);
00548     N_free_array_3d(data->disp_yy);
00549     N_free_array_3d(data->disp_zz);
00550     N_free_array_3d(data->disp_xy);
00551     N_free_array_3d(data->disp_xz);
00552     N_free_array_3d(data->disp_yz);
00553 
00554     G_free(data);
00555 
00556     data = NULL;
00557 
00558     return;
00559 }
00560 
00561 /* ************************************************************************* *
00562  * ************************************************************************* *
00563  * ************************************************************************* */
00570 void N_free_solute_transport_data2d(N_solute_transport_data2d * data)
00571 {
00572     N_free_array_2d(data->c);
00573     N_free_array_2d(data->c_start);
00574     N_free_array_2d(data->status);
00575     N_free_array_2d(data->diff_x);
00576     N_free_array_2d(data->diff_y);
00577     N_free_array_2d(data->q);
00578     N_free_array_2d(data->cs);
00579     N_free_array_2d(data->R);
00580     N_free_array_2d(data->nf);
00581     N_free_array_2d(data->cin);
00582     N_free_array_2d(data->top);
00583     N_free_array_2d(data->bottom);
00584 
00585     N_free_array_2d(data->disp_xx);
00586     N_free_array_2d(data->disp_yy);
00587     N_free_array_2d(data->disp_xy);
00588 
00589     G_free(data);
00590 
00591     data = NULL;
00592 
00593     return;
00594 }
00595 
00612 void N_calc_solute_transport_transmission_2d(N_solute_transport_data2d * data)
00613 {
00614     int i, j, count = 1;
00615 
00616     int cols, rows;
00617 
00618     double c;
00619 
00620     N_gradient_2d grad;
00621 
00622     cols = data->grad->cols;
00623     rows = data->grad->rows;
00624 
00625     G_debug(2,
00626             "N_calc_solute_transport_transmission_2d: calculating transmission boundary");
00627 
00628     for (j = 0; j < rows; j++) {
00629         for (i = 0; i < cols; i++) {
00630             if (N_get_array_2d_d_value(data->status, i, j) ==
00631                 N_CELL_TRANSMISSION) {
00632                 count = 0;
00633                 /*get the gradient neighbours */
00634                 N_get_gradient_2d(data->grad, &grad, i, j);
00635                 c = 0;
00636                 /*
00637                    c = N_get_array_2d_d_value(data->c_start, i, j);
00638                    if(c > 0)
00639                    count++;
00640                  */
00641 
00642                 if (grad.WC > 0 &&
00643                     !N_is_array_2d_value_null(data->c, i - 1, j)) {
00644                     c += N_get_array_2d_d_value(data->c, i - 1, j);
00645                     count++;
00646                 }
00647                 if (grad.EC < 0 &&
00648                     !N_is_array_2d_value_null(data->c, i + 1, j)) {
00649                     c += N_get_array_2d_d_value(data->c, i + 1, j);
00650                     count++;
00651                 }
00652                 if (grad.NC < 0 &&
00653                     !N_is_array_2d_value_null(data->c, i, j - 1)) {
00654                     c += N_get_array_2d_d_value(data->c, i, j - 1);
00655                     count++;
00656                 }
00657                 if (grad.SC > 0 &&
00658                     !N_is_array_2d_value_null(data->c, i, j + 1)) {
00659                     c += N_get_array_2d_d_value(data->c, i, j + 1);
00660                     count++;
00661                 }
00662                 if (count != 0)
00663                     c = c / (double)count;
00664                 /*make sure it is not NAN */
00665                 if (c > 0 || c == 0 || c < 0)
00666                     N_put_array_2d_d_value(data->c_start, i, j, c);
00667             }
00668         }
00669     }
00670 
00671     return;
00672 }
00673 
00688 void N_calc_solute_transport_disptensor_2d(N_solute_transport_data2d * data)
00689 {
00690     int i, j;
00691 
00692     int cols, rows;
00693 
00694     double vx, vy, vv;
00695 
00696     double disp_xx, disp_yy, disp_xy;
00697 
00698     N_gradient_2d grad;
00699 
00700     cols = data->grad->cols;
00701     rows = data->grad->rows;
00702 
00703     G_debug(2,
00704             "N_calc_solute_transport_disptensor_2d: calculating the dispersivity tensor");
00705 
00706     for (j = 0; j < rows; j++) {
00707         for (i = 0; i < cols; i++) {
00708 
00709             disp_xx = 0;
00710             disp_yy = 0;
00711             disp_xy = 0;
00712 
00713             /*get the gradient neighbours */
00714             N_get_gradient_2d(data->grad, &grad, i, j);
00715             vx = (grad.WC + grad.EC) / 2;
00716             vy = (grad.NC + grad.SC) / 2;
00717             vv = sqrt(vx * vx + vy * vy);
00718 
00719             if (vv != 0) {
00720                 disp_xx = data->al * vx * vx / vv + data->at * vy * vy / vv;
00721                 disp_yy = data->at * vx * vx / vv + data->al * vy * vy / vv;
00722                 disp_xy = (data->al - data->at) * vx * vy / vv;
00723             }
00724 
00725             G_debug(5,
00726                     "N_calc_solute_transport_disptensor_2d: [%i][%i] disp_xx %g disp_yy %g disp_xy %g",
00727                     i, j, disp_xx, disp_yy, disp_xy);
00728             N_put_array_2d_d_value(data->disp_xx, i, j, disp_xx);
00729             N_put_array_2d_d_value(data->disp_yy, i, j, disp_yy);
00730             N_put_array_2d_d_value(data->disp_xy, i, j, disp_xy);
00731         }
00732     }
00733 
00734     return;
00735 }
00736 
00751 void N_calc_solute_transport_disptensor_3d(N_solute_transport_data3d * data)
00752 {
00753     int i, j, k;
00754 
00755     int cols, rows, depths;
00756 
00757     double vx, vy, vz, vv;
00758 
00759     double disp_xx, disp_yy, disp_zz, disp_xy, disp_xz, disp_yz;
00760 
00761     N_gradient_3d grad;
00762 
00763     cols = data->grad->cols;
00764     rows = data->grad->rows;
00765     depths = data->grad->depths;
00766 
00767     G_debug(2,
00768             "N_calc_solute_transport_disptensor_3d: calculating the dispersivity tensor");
00769 
00770     for (k = 0; k < depths; k++) {
00771         for (j = 0; j < rows; j++) {
00772             for (i = 0; i < cols; i++) {
00773                 disp_xx = 0;
00774                 disp_yy = 0;
00775                 disp_zz = 0;
00776                 disp_xy = 0;
00777                 disp_xz = 0;
00778                 disp_yz = 0;
00779 
00780                 /*get the gradient neighbours */
00781                 N_get_gradient_3d(data->grad, &grad, i, j, k);
00782                 vx = (grad.WC + grad.EC) / 2;
00783                 vy = (grad.NC + grad.SC) / 2;
00784                 vz = (grad.BC + grad.TC) / 2;
00785                 vv = sqrt(vx * vx + vy * vy + vz * vz);
00786 
00787                 if (vv != 0) {
00788                     disp_xx =
00789                         data->al * vx * vx / vv + data->at * vy * vy / vv +
00790                         data->at * vz * vz / vv;
00791                     disp_yy =
00792                         data->at * vx * vx / vv + data->al * vy * vy / vv +
00793                         data->at * vz * vz / vv;
00794                     disp_zz =
00795                         data->at * vx * vx / vv + data->at * vy * vy / vv +
00796                         data->al * vz * vz / vv;
00797                     disp_xy = (data->al - data->at) * vx * vy / vv;
00798                     disp_xz = (data->al - data->at) * vx * vz / vv;
00799                     disp_yz = (data->al - data->at) * vy * vz / vv;
00800                 }
00801 
00802                 G_debug(5,
00803                         "N_calc_solute_transport_disptensor_3d: [%i][%i][%i] disp_xx %g disp_yy %g disp_zz %g  disp_xy %g disp_xz %g disp_yz %g ",
00804                         i, j, k, disp_xx, disp_yy, disp_zz, disp_xy, disp_xz,
00805                         disp_yz);
00806                 N_put_array_3d_d_value(data->disp_xx, i, j, k, disp_xx);
00807                 N_put_array_3d_d_value(data->disp_yy, i, j, k, disp_yy);
00808                 N_put_array_3d_d_value(data->disp_zz, i, j, k, disp_zz);
00809                 N_put_array_3d_d_value(data->disp_xy, i, j, k, disp_xy);
00810                 N_put_array_3d_d_value(data->disp_xz, i, j, k, disp_xz);
00811                 N_put_array_3d_d_value(data->disp_yz, i, j, k, disp_yz);
00812             }
00813         }
00814     }
00815 
00816     return;
00817 }