TheBoussinesqModel  3.2.1
 All Data Structures Files Functions Variables Typedefs Macros Pages
b_solver.c File Reference
#include "turtle.h"
#include "t_alloc.h"
#include "t_datamanipulation.h"
#include "t_utilities.h"
#include "rw_maps.h"
#include "doublevector_utilities.h"
#include "preconditioned_conjugate_gradient.h"
#include "keywords_file_b.h"
#include "geometry.h"
#include "g_raster2plvector.h"
#include "bigcells2.h"
#include "geometry2.h"
#include "b_utilities.h"
#include "b_volumes.h"
#include "b_sources.h"
#include "b_v_advection.h"
#include "b_solver.h"

Go to the source code of this file.

Macros

#define MIN_T_DIAG_NEGATIVE_INDEX   -20
 
#define MIN_T_DIAG   1.0e-10
 

Functions

double t_st_operator_element (long i, DOUBLEVECTOR *eta)
 
double t_st_operator_element_no_dirichlet (long i, DOUBLEVECTOR *eta)
 
double t_st_operator_element_subs (long i, DOUBLEVECTOR *eta, double cond_dirichlet)
 
double b_smatrix_element (long i, DOUBLEVECTOR *x)
 
long Newton_convergence (DOUBLEVECTOR *x_temp, DOUBLEVECTOR *be, DOUBLEVECTOR *be0)
 
int time_loop (short print, int(*write_output)(void *v1, void *v2))
 
double water_surface_elevation_mean (double eta1, double eta2)
 
int write_map_results (void *output, void *SSSS)
 
OUTPUT_FILENAMESnew_output_filenames (short print)
 
void free_OUTPUT_FILENAMES (OUTPUT_FILENAMES *fn)
 

Variables

STRINGBINfilenames
 
DOUBLESQUARE_GRIDdsq
 
DOUBLE_GRIDdgrid
 
DOUBLERASTER_MAPdraster
 
char * wpath
 
PARAMparam
 
FLAGflag
 
DOUBLEVECTORelevation_bottom_fine
 
DOUBLEVECTORelevation_bottom_coarse
 
DOUBLEVECTORelevation_bottom_flines
 
DOUBLEVECTORelevation_bottom_bottom
 
DOUBLEVECTORporosity_fine
 
DOUBLEVECTORoutlet_coefficient
 
DOUBLEVECTORwater_surface_elevation
 
DOUBLEVECTORwater_mass_error
 
DOUBLEVECTORwater_depth
 
DOUBLEVECTORsurface_water_velocity
 
DOUBLEVECTORF1_wet_vert_area
 
DOUBLEVECTOReta_v
 
DOUBLEVECTORsources
 
DOUBLEVECTORdirichlet
 
DOUBLEVECTORt_diagonal
 
DOUBLEVECTORt_diagonal_no_dirichlet
 
S_TIMESs_times
 
S_TIMESdirichlet_times
 

Detailed Description

Author
Emanuele Cordano
Attention

This file is part of Boussinesq.

Boussinesq is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

Boussinesq is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with Boussinesq. If not, see http://www.gnu.org/licenses/.

Definition in file b_solver.c.

Macro Definition Documentation

#define MIN_T_DIAG   1.0e-10

Definition at line 50 of file b_solver.c.

Referenced by Newton_convergence().

#define MIN_T_DIAG_NEGATIVE_INDEX   -20

Definition at line 48 of file b_solver.c.

Referenced by b_smatrix_element(), and Newton_convergence().

Function Documentation

double b_smatrix_element ( long  i,
DOUBLEVECTOR x 
)
Author
Emanuele Cordano
Date
29 June 2009
Parameters
i- index of the polygon at which t_st_operator_element is calculated
eta- (DOUBLEVECTOR *) vector of the unknowns

Function inverted and solved with a preconditioned conjugate gradient method obtained as follows:

\[ \left[ \mathbf{T} \cdot \eta^{n+1} \right]_{i}= -\Delta t^{n} \, \sum_{j \in S_i} \, K_S \, A_j(\eta^n_{i,m(i,j)})\frac{x_{m(i,j)}^{n+1}-x_{i}^{n+1}}{\delta_{i,m(i,j)}}+ Ws(^m\eta^{n+1})x_i \]

which is the t_st_operator_element() function plus x_i multiplied the wet area of the i-th pixel. This function has the same signature of t_st_operator_element() Please see the documentation in pdf with the algebraic passages for further details.

Definition at line 274 of file b_solver.c.

References PARAM::deta, MIN_T_DIAG_NEGATIVE_INDEX, PARAM::null_dirichlet, t_st_operator_element_no_dirichlet(), and wet_area().

Referenced by Newton_convergence().

void free_OUTPUT_FILENAMES ( OUTPUT_FILENAMES fn)

Definition at line 801 of file b_solver.c.

References OUTPUT_FILENAMES::file_error, and OUTPUT_FILENAMES::file_result.

Referenced by time_loop().

long Newton_convergence ( DOUBLEVECTOR x_temp,
DOUBLEVECTOR be,
DOUBLEVECTOR be0 
)
Author
Emanuele Cordano
Date
20 April 2009
Parameters
x_term(DOUBLEVECTOR *) x_term - doublevector of residuals (unkowns)
be(DOUBLEVECTOR *) be - known term of the equation in Newton Itaration Method (balance equation);
be0(DOUBLEVECTOR *) be - known term of the equation in .... (balance equation);
  • This functions contains the iterative the Newton-like iterative according to Casulli, 2008 :

    \[ ^{m+1}\eta^{n+1}= ^{m}\eta^{n+1}-\left[\mathbf{Ws}( ^{m}\eta^{n+1})+\mathbf{T}\right]^{-1}\left[\mathbf{V}( ^{m}\eta^{n+1})+\mathbf{T}\cdot^{m}\eta^{n+1} -\mathbf{b3}\right] \]

    where $m$ is the reiteration level. Please see the pdf documentation of the analytical passages for the documentation.

Definition at line 367 of file b_solver.c.

References b_smatrix_element(), free_doublevector(), jacobi_preconditioned_conjugate_gradient_search(), join_strings(), max_doublevector(), PARAM::max_errors, MIN_T_DIAG, MIN_T_DIAG_NEGATIVE_INDEX, new_doublevector(), DOUBLEVECTOR::nh, DOUBLEVECTOR::nl, PARAM::null_dirichlet, t_st_operator_element(), volume(), wpath, and PARAM::x_temp_adm.

Referenced by time_loop().

double t_st_operator_element ( long  i,
DOUBLEVECTOR eta 
)
Author
Emanuele Cordano
Date
October 2009

Definition at line 75 of file b_solver.c.

References MAX_ELEVATION_VALUE, t_st_advection_operator_element(), and t_st_operator_element_subs().

Referenced by Newton_convergence(), and time_loop().

double t_st_operator_element_no_dirichlet ( long  i,
DOUBLEVECTOR eta 
)
Author
Emanuele Cordano
Date
May 2010

Definition at line 87 of file b_solver.c.

References PARAM::null_dirichlet, t_st_advection_operator_element(), and t_st_operator_element_subs().

Referenced by b_smatrix_element(), and time_loop().

double t_st_operator_element_subs ( long  i,
DOUBLEVECTOR eta,
double  cond_dirichlet 
)
Author
Emanuele Cordano
Date
June 2009
Parameters
i- index of the polygon at which t_st_operator_element is calculated
eta- (DOUBLEVECTOR *) vector of the unknowns
cond_dirichlet- (double) (=MAX_ELEVATION_VALUE) in case dirichlet value are neglected or (=param->null_dirichlet) othrwise (added by ec on 20100517)
       This function calculates the i-th element of the following vector
       where eta is the vector of the independent variables

\[ \left[ \mathbf{T} \cdot \eta^{n+1} \right]_{i}= -\Delta t^{n} \, \sum_{j \in S_i} \, K_S \, A_j(\eta^n_{i,m(i,j)})\frac{\eta_{m(i,j)}^{n+1}-\eta_{i}^{n+1}}{\delta_{i,m(i,j)}} \]

Please see the pdf documentation with the algebraic passages for major details.

Definition at line 99 of file b_solver.c.

References GRID::boundary_indicator, DOUBLE_GRID::coarse, polygon_connection_attributes::connections, polygon_connection_attributes::d_connections, PARAM::dt, POLYGON::edge_indices, POLYGONVECTOR::element, polygon_connection_attribute_array::element, PARAM::Ks, GRID::links, LONGVECTOR::nl, GRID::polygons, vertical_area_subs(), and water_surface_elevation_mean().

Referenced by t_st_operator_element(), and t_st_operator_element_no_dirichlet().

int time_loop ( short  print,
int(*)(void *v1, void *v2)  write_output 
)
Author
Emanuele Cordano
Date
22 April 2009
Parameters
(short)- print FLuidTurtle Parametesrs
int(*write_output)(void *v1, void *v2) - function which writes the output for every requested time steps
Returns
value 0 in case of success, -1 otherwise

In this routin the discretized system equation:

\[ \mathbf{V}(\eta^{n+1})+\mathbf{T} \cdot \eta^{n+1} = \mathbf{b3} \]

is solved. Please see the meanings of the symbols in the pdf documentation with algebraic passages. The nonlinear equation system is solved with a Newtton-like

Definition at line 548 of file b_solver.c.

References POLYGON::area2D, b_advection(), DOUBLE_GRID::coarse, PARAM::dt, PARAM::dt_print, POLYGONVECTOR::element, free_doublevector(), free_OUTPUT_FILENAMES(), free_s_times(), get_diagonal(), get_dirichletsnode(), get_s_times(), get_sources(), I_DIRICHLETTIMES_FT_FILE, I_TIMES_FT_FILE, new_doublevector(), new_output_filenames(), Newton_convergence(), DOUBLEVECTOR::nh, DOUBLEVECTOR::nl, PARAM::null_dirichlet, GRID::polygons, q_discharge_from_outlet_cell(), PARAM::t, PARAM::t_end, t_st_operator_element(), t_st_operator_element_no_dirichlet(), PARAM::t_start, update_F1_wet_vert_area(), update_velocity(), volume(), and write_suffix().

Referenced by main().

double water_surface_elevation_mean ( double  eta1,
double  eta2 
)
Author
Emanuele Cordano
Date
June 2010

Definition at line 714 of file b_solver.c.

References FLAG::arithmetic_mean0.

Referenced by b_advection(), t_st_advection_operator_element(), t_st_operator_element_subs(), update_F1_wet_vert_area(), and update_velocity().

Variable Documentation

DOUBLE_GRID* dgrid

Definition at line 126 of file main.c.

Referenced by new_double_grid_from_doublesquare_grid().

DOUBLEVECTOR * dirichlet

coefficient F1 for surface flow as reported in Casulli,2008

Definition at line 71 of file b_solver.c.

S_TIMES * dirichlet_times

Definition at line 72 of file b_solver.c.

DOUBLERASTER_MAP* draster

Definition at line 128 of file main.c.

Referenced by get_doubleraster_map().

Definition at line 125 of file main.c.

Referenced by get_doublesquare_grid(), and read_doublesquare_grid().

DOUBLEVECTOR* elevation_bottom_bottom

map of the bottom elevation defined in the lines of the fine grid

Definition at line 134 of file main.c.

DOUBLEVECTOR* elevation_bottom_coarse

map of the bottom elevation (fine grid)

Definition at line 132 of file main.c.

DOUBLEVECTOR* elevation_bottom_fine

Definition at line 131 of file main.c.

DOUBLEVECTOR* elevation_bottom_flines

map of the bottom elevation (coarse grid)

Definition at line 133 of file main.c.

DOUBLEVECTOR* eta_v

coefficient F1 for surface flow as reported in Casulli,2008

Definition at line 71 of file b_solver.c.

DOUBLEVECTOR* F1_wet_vert_area

map of the surface velocity defined in the lines (velocity is positive when going kp1 to kp2 polygons with kp1>kp2 )

map of the (wet area) surface velocity defined in the lines (velocity is positive when going kp1 to kp2 polygons with kp1>kp2 )

Definition at line 147 of file main.c.

STRINGBIN* filenames

Definition at line 124 of file main.c.

FLAG* flag

Definition at line 149 of file main.c.

DOUBLEVECTOR* outlet_coefficient

map of porosity defined on the pixels of the fine grid map of the coefficient of the rating curve in the outlet q_discharge=C*h_sup^m

map of porosity defined on the pixels of the fine grid

Definition at line 138 of file main.c.

PARAM* param

coefficient F1 for surface flow as reported in Casulli,2008

Definition at line 148 of file main.c.

DOUBLEVECTOR* porosity_fine

map of the elvetion of the bottom of he coarse cells

map of the surface terrain elevation (fine grid lines)

Definition at line 137 of file main.c.

S_TIMES* s_times

eta_v and sources and dirichlet are defined as a global value, t_diagonal is the diagonal of t_st_element tensor

Definition at line 72 of file b_solver.c.

DOUBLEVECTOR * sources

Definition at line 71 of file b_solver.c.

DOUBLEVECTOR* surface_water_velocity

Definition at line 146 of file main.c.

DOUBLEVECTOR * t_diagonal

Definition at line 71 of file b_solver.c.

DOUBLEVECTOR * t_diagonal_no_dirichlet

Definition at line 71 of file b_solver.c.

DOUBLEVECTOR* water_depth

map of instantaneous water mass error

Definition at line 144 of file main.c.

DOUBLEVECTOR* water_mass_error

map of instantaneous water surface on the pixels of the coarse grid

Definition at line 143 of file main.c.

DOUBLEVECTOR* water_surface_elevation

Definition at line 142 of file main.c.

char* wpath

Definition at line 151 of file main.c.

Referenced by main(), and Newton_convergence().