cfs_tpfa_residual.h
Go to the documentation of this file.
1 /*
2  Copyright 2010 SINTEF ICT, Applied Mathematics.
3 
4  This file is part of the Open Porous Media project (OPM).
5 
6  OPM is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  OPM is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with OPM. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #ifndef OPM_CFS_TPFA_HEADER_INCLUDED
21 #define OPM_CFS_TPFA_HEADER_INCLUDED
22 
23 #include <opm/core/grid.h>
24 #include <opm/core/wells.h>
25 
27 
59 #ifdef __cplusplus
60 extern "C" {
61 #endif
62 
63 struct cfs_tpfa_res_impl;
64 struct CSRMatrix;
66 
76  struct Wells *W ;
77 
82 };
83 
84 
90  struct compr_src *src ;
91 };
92 
93 
99  struct CSRMatrix *J;
100  double *F;
102  struct cfs_tpfa_res_impl *pimpl;
103 };
104 
105 
123 struct cfs_tpfa_res_data *
125  struct cfs_tpfa_res_wells *wells ,
126  int nphases);
127 
128 
142 void
144 
145 
182 int
184  double dt,
185  struct cfs_tpfa_res_forces *forces,
186  const double *zc,
187  struct compr_quantities_gen *cq,
188  const double *trans,
189  const double *gravcap_f,
190  const double *cpress,
191  const double *wpress,
192  const double *porevol,
193  struct cfs_tpfa_res_data *h);
194 
195 
237 int
239  struct UnstructuredGrid *G,
240  double dt,
241  struct cfs_tpfa_res_forces *forces,
242  const double *zc,
243  struct compr_quantities_gen *cq,
244  const double *trans,
245  const double *gravcap_f,
246  const double *cpress,
247  const double *wpress,
248  const double *porevol,
249  const double *porevol0,
250  const double *rock_comp,
251  struct cfs_tpfa_res_data *h);
252 
253 
285 void
287  struct cfs_tpfa_res_forces *forces ,
288  int np ,
289  const double *trans ,
290  const double *pmobc ,
291  const double *pmobf ,
292  const double *gravcap_f,
293  const double *cpress ,
294  const double *wpress ,
295  double *fflux ,
296  double *wflux );
297 
298 
323 void
325  int np,
326  const double *htrans,
327  const double *pmobf,
328  const double *gravcap_f,
329  struct cfs_tpfa_res_data *h,
330  const double *cpress,
331  const double *fflux,
332  double *fpress);
333 
334 #if 0
335 void
337  int np,
338  struct cfs_tpfa_data *h,
339  double *masstrans_f);
340 
341 void
343  int np,
344  struct cfs_tpfa_data *h,
345  double *gravtrans_f);
346 
347 double
349  struct compr_quantities *cq,
350  const double *trans,
351  const double *porevol,
352  struct cfs_tpfa_data *h,
353  const double *dpmobf,
354  const double *surf_dens,
355  const double *gravity);
356 
357 void
359  well_t *W,
360  struct completion_data *wdata,
361  int np,
362  double dt,
363  const double *porevol,
364  struct cfs_tpfa_data *h,
365  double *surf_vol);
366 
367 #endif
368 
369 #ifdef __cplusplus
370 }
371 #endif
372 
373 #endif /* OPM_CFS_TPFA_HEADER_INCLUDED */
Definition: compr_source.h:52
struct cfs_tpfa_res_data * cfs_tpfa_res_construct(struct UnstructuredGrid *G, struct cfs_tpfa_res_wells *wells, int nphases)
Definition: wells.h:50
int cfs_tpfa_res_assemble(struct UnstructuredGrid *G, double dt, struct cfs_tpfa_res_forces *forces, const double *zc, struct compr_quantities_gen *cq, const double *trans, const double *gravcap_f, const double *cpress, const double *wpress, const double *porevol, struct cfs_tpfa_res_data *h)
void cfs_tpfa_retrieve_masstrans(struct UnstructuredGrid *G, int np, struct cfs_tpfa_data *h, double *masstrans_f)
struct Wells * W
Definition: cfs_tpfa_residual.h:76
Definition: sparse_sys.h:38
Definition: grid.h:98
struct cfs_tpfa_res_impl * pimpl
Definition: cfs_tpfa_residual.h:102
Definition: compr_quant.h:31
struct compr_src * src
Definition: cfs_tpfa_residual.h:90
Definition: cfs_tpfa_residual.h:98
Definition: cfs_tpfa_residual.h:71
int cfs_tpfa_res_comprock_assemble(struct UnstructuredGrid *G, double dt, struct cfs_tpfa_res_forces *forces, const double *zc, struct compr_quantities_gen *cq, const double *trans, const double *gravcap_f, const double *cpress, const double *wpress, const double *porevol, const double *porevol0, const double *rock_comp, struct cfs_tpfa_res_data *h)
Definition: cfs_tpfa.h:35
Definition: wells.h:132
Definition: compr_quant_general.h:40
double cfs_tpfa_impes_maxtime(struct UnstructuredGrid *G, struct compr_quantities *cq, const double *trans, const double *porevol, struct cfs_tpfa_data *h, const double *dpmobf, const double *surf_dens, const double *gravity)
void cfs_tpfa_res_flux(struct UnstructuredGrid *G, struct cfs_tpfa_res_forces *forces, int np, const double *trans, const double *pmobc, const double *pmobf, const double *gravcap_f, const double *cpress, const double *wpress, double *fflux, double *wflux)
void cfs_tpfa_retrieve_gravtrans(struct UnstructuredGrid *G, int np, struct cfs_tpfa_data *h, double *gravtrans_f)
struct cfs_tpfa_res_wells * wells
Definition: cfs_tpfa_residual.h:89
struct CompletionData * data
Definition: cfs_tpfa_residual.h:81
const double gravity
Definition: Units.hpp:120
struct CSRMatrix * J
Definition: cfs_tpfa_residual.h:99
void cfs_tpfa_res_fpress(struct UnstructuredGrid *G, int np, const double *htrans, const double *pmobf, const double *gravcap_f, struct cfs_tpfa_res_data *h, const double *cpress, const double *fflux, double *fpress)
double * F
Definition: cfs_tpfa_residual.h:100
Definition: legacy_well.h:71
void cfs_tpfa_res_destroy(struct cfs_tpfa_res_data *h)
Definition: legacy_well.h:52
Definition: cfs_tpfa_residual.h:88
void cfs_tpfa_expl_mass_transport(struct UnstructuredGrid *G, well_t *W, struct completion_data *wdata, int np, double dt, const double *porevol, struct cfs_tpfa_data *h, double *surf_vol)