TrioCFD 1.9.8
TrioCFD documentation
Loading...
Searching...
No Matches
Operateur_IJK_base.cpp
1/****************************************************************************
2* Copyright (c) 2015 - 2016, CEA
3* All rights reserved.
4*
5* Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
6* 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
7* 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution.
8* 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission.
9*
10* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
11* IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
12* OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
13*
14*****************************************************************************/
15
16#include <Operateur_IJK_base.h>
17#include <IJK_Field_vector.h>
18#include <Cut_cell_FT_Disc.h>
19// Definition of the fluxes
20//
21// control volumes for velocity in x direction:
22// this face_i(i,j,k) is located at left of element(i,j,k)
23//
24// -------------------------
25// . | .
26// . | .
27// . | .
28// . | .
29// . |-----> .
30// . | .
31// . | .
32// . | .
33// -------------------------
34// face(i,j,k) receives fluxes from all neighbour faces:
35// flux_i(i,j,k) at left
36// flux_j(i,j,k)
37// flux_k(i,j,k)
38// flux_i(i+1,j,k) at right
39// ...
40
41static inline void swap_data(IJK_Field_local_double*& a, IJK_Field_local_double*& b)
42{
43 IJK_Field_local_double *tmp=a;
44 a=b;
45 b=tmp;
46}
47
48/*
49 * Operateur_IJK_faces_base
50 */
51
52Implemente_base(Operateur_IJK_faces_base_double,"Operateur_IJK_faces_base_double", Objet_U);
53
55{
56 // Objet_U::printOn(os);
57 return os;
58}
59
61{
62 // Objet_U::readOn(is);
63 return is;
64}
65
66// Compute the maximum stable timestep for a convection scheme on the local processor.
67// The global stability timestep is the minimum of all values on all processors.
68// The ghost velocities must be uptodate on the right each subdomain.
69double Operateur_IJK_faces_base_double::compute_dtstab_convection_local(IJK_Field_double& vx, IJK_Field_double& vy, IJK_Field_double& vz)
70{
71 double inverse_dt_conv = 0.;
72 const Domaine_IJK& geom = vx.get_domaine();
73 const int ni = geom.get_nb_elem_local(DIRECTION_I);
74 const int nj = geom.get_nb_elem_local(DIRECTION_J);
75 const int nk = geom.get_nb_elem_local(DIRECTION_K);
76 const double dx = geom.get_constant_delta(DIRECTION_I);
77 const double dy = geom.get_constant_delta(DIRECTION_J);
78 const ArrOfDouble& delta_z = geom.get_delta(DIRECTION_K);
79 const int k_offset = geom.get_offset_local(DIRECTION_K);
80
81 for (int k = 0 ; k < nk ; k++)
82 {
83 // pointeur pour optimisation du calcul des voisins ( marche pour les derivees)
84 ConstIJK_double_ptr ptr_vx(vx, 0, 0, k);
85 ConstIJK_double_ptr ptr_vy(vy, 0, 0, k);
86 ConstIJK_double_ptr ptr_vz(vz, 0, 0, k);
87 const double dz = delta_z[k + k_offset];
88 const double inverse_dx = 1. / dx;
89 const double inverse_dy = 1. / dy;
90 const double inverse_dz = 1. / dz;
91
92 for (int j = 0 ; j < nj ; j++)
93 {
94 for (int i = 0 ; i < ni ; i++)
95 {
96 // pas de temps de convection
97 double vx0 , vx1 , vy0 , vy1 , vz0 , vz1 ; // vitesses face gauche en 0 et face droite en 1
98 ptr_vx.get_center_right(DIRECTION::X, i, vx0, vx1);
99 ptr_vy.get_center_right(DIRECTION::Y, i, vy0, vy1);
100 ptr_vz.get_center_right(DIRECTION::Z, i, vz0, vz1);
101 double inverse_dt_conv_loc =
102 (std::max(vx0, 0.) - std::min(vx1, 0.)) * inverse_dx
103 + (std::max(vy0, 0.) - std::min(vy1, 0.)) * inverse_dy
104 + (std::max(vz0, 0.) - std::min(vz1, 0.)) * inverse_dz;
105
106 inverse_dt_conv = std::max(inverse_dt_conv, inverse_dt_conv_loc);
107 }
108 if ( j < nj-1 )
109 {
110 ptr_vx.next_j();
111 ptr_vy.next_j();
112 ptr_vz.next_j();
113 }
114 }
115 }
116 inverse_dt_conv = std::max(1e-20, inverse_dt_conv);
117 return 1. / inverse_dt_conv;
118}
119
120void Operateur_IJK_faces_base_double::compute_set(IJK_Field_double& dvx, IJK_Field_double& dvy, IJK_Field_double& dvz)
121{
122 compute_(dvx, dvy, dvz, 0 /* setting */);
123}
124
125void Operateur_IJK_faces_base_double::compute_add(IJK_Field_double& dvx, IJK_Field_double& dvy, IJK_Field_double& dvz)
126{
127 compute_(dvx, dvy, dvz, 1 /* adding */);
128}
129
130void Operateur_IJK_faces_base_double::compute_(IJK_Field_double& dvx, IJK_Field_double& dvy, IJK_Field_double& dvz, bool add)
131{
132 IJK_Field_local_double storage;
133 storage.allocate(dvx.ni()+1, dvy.nj()+1, 6, 0);
134 IJK_Field_local_double tmp[6];
135 for (int i = 0; i < 6; i++)
136 tmp[i].ref_ij(storage, i);
137 IJK_Field_local_double *flux_vx_zmin = &tmp[0]; // flux pour la composante x a travers la face zmin des volumes de controles
138 IJK_Field_local_double *flux_vy_zmin = &tmp[1]; // flux pour la composante y
139 IJK_Field_local_double *flux_vz_zmin = &tmp[2]; // flux pour la composante z
140 IJK_Field_local_double *flux_zmax = &tmp[3]; // flux a travers la face zmax des volumes de controles
141 IJK_Field_local_double * const flux_x = &tmp[4]; // pointer will not change, content will change
142 IJK_Field_local_double * const flux_y = &tmp[5];
143
144 compute_flux_z_vx(*flux_vx_zmin, 0);
145 compute_flux_z_vy(*flux_vy_zmin, 0);
146 compute_flux_z_vz(*flux_vz_zmin, 0);
147
148 correct_flux(flux_vx_zmin, 0, 2, 0);
149 correct_flux(flux_vy_zmin, 0, 2, 1);
150 correct_flux(flux_vz_zmin, 0, 2, 2);
151
152 // At the bottom end of the domain, there is more layer of unknowns for vz:
153 const int kmax = std::max(std::max(dvx.nk(), dvy.nk()), dvz.nk());
154 for (int k = 0; k < kmax ; k++)
155 {
156 if (k < dvz.nk())
157 {
158
159 compute_flux_x_vz(*flux_x, k);
160 compute_flux_y_vz(*flux_y, k);
161 compute_flux_z_vz(*flux_zmax, k+1);
162
163 correct_flux(flux_x, k, 0, 2);
164 correct_flux(flux_y, k, 1, 2);
165 correct_flux(flux_zmax, k+1, 2, 2);
166
167 Operator_IJK_div(*flux_x, *flux_y, *flux_vz_zmin, *flux_zmax, dvz, k, add);
168
170 swap_data(flux_vz_zmin, flux_zmax); // conserve le flux en z pour la couche z suivante
171 }
172 if (k < dvx.nk())
173 {
174 compute_flux_x_vx(*flux_x, k);
175 compute_flux_y_vx(*flux_y, k);
176 compute_flux_z_vx(*flux_zmax, k+1);
177
178 correct_flux(flux_x, k, 0, 0);
179 correct_flux(flux_y, k, 1, 0);
180 correct_flux(flux_zmax, k+1, 2, 0);
181
182 Operator_IJK_div(*flux_x, *flux_y, *flux_vx_zmin, *flux_zmax, dvx, k, add);
183
185 swap_data(flux_vx_zmin, flux_zmax); // conserve le flux en z pour la couche z suivante
186 }
187 if (k < dvy.nk())
188 {
189 compute_flux_x_vy(*flux_x, k);
190 compute_flux_y_vy(*flux_y, k);
191 compute_flux_z_vy(*flux_zmax, k+1);
192
193 correct_flux(flux_x, k, 0, 1);
194 correct_flux(flux_y, k, 1, 1);
195 correct_flux(flux_zmax, k+1, 2, 1);
196
197 Operator_IJK_div(*flux_x, *flux_y, *flux_vy_zmin, *flux_zmax, dvy, k, add);
198
200 swap_data(flux_vy_zmin, flux_zmax); // conserve le flux en z pour la couche z suivante
201 }
202 }
203}
204
205// Compute the "- divergence" of the given flux, integrated over control volume
206// We assume than flux contain the integral of the flux on the faces of the control volumes
207// For "div_set" method:
208// resu(i,j,k) = flux_x(i,j) - flux_x(i+1,j) + flux_y(i,j) - flux_y(i,j+1) + flux_zmin(i,j) - flux_zmax(i,j)
209// For "div_add" method: same but
210// resu(i,j,k) += ...
211void Operateur_IJK_faces_base_double::Operator_IJK_div(const IJK_Field_local_double& flux_x, const IJK_Field_local_double& flux_y,
212 const IJK_Field_local_double& flux_zmin, const IJK_Field_local_double& flux_zmax,
213 IJK_Field_double& resu, int k_layer, bool add)
214{
215 ConstIJK_double_ptr fx(flux_x, 0, 0, 0);
216 ConstIJK_double_ptr fy(flux_y, 0, 0, 0);
217 ConstIJK_double_ptr fzmin(flux_zmin, 0, 0, 0);
218 ConstIJK_double_ptr fzmax(flux_zmax, 0, 0, 0);
219 IJK_double_ptr resu_ptr(resu, 0, 0, k_layer);
220
221
222 const int imax = resu.ni();
223 const int jmax = resu.nj();
224 const int vsize = Simd_double::size();
225 for (int j = 0; ; j++)
226 {
227 for (int i = 0; i < imax; i += vsize)
228 {
229 Simd_double r, a, b;
230 fx.get_center_right(DIRECTION::X, i, a, b);
231 r = a-b;
232
233 fy.get_center_right(DIRECTION::Y,i, a, b);
234 r += a-b;
235
236 fzmin.get_center(i, a);
237 fzmax.get_center(i, b);
238 r += a-b;
239
240 if (add)
241 {
242 resu_ptr.get_center(i, a);
243 r += a;
244 }
245 resu_ptr.put_val(i, r);
246 }
247 // do not execute end_iloop at last iteration (because of assert on valid j+1)
248 if (j+1==jmax)
249 break;
250 fx.next_j();
251 fy.next_j();
252 fzmin.next_j();
253 fzmax.next_j();
254 resu_ptr.next_j();
255 }
256
257}
258
259void Operateur_IJK_faces_base_double::correct_flux(IJK_Field_local_double *const flux, const int k_layer, const int dir, const int compo)
260{
261}
262
263
264/*
265 * Operateur_IJK_elem_base
266 */
267
268Implemente_base(Operateur_IJK_elem_base_double,"Operateur_IJK_elem_base_double",Objet_U);
269
271{
272 // Objet_U::printOn(os);
273 return os;
274}
275
277{
278 // Objet_U::readOn(is);
279 return is;
280}
281
282void Operateur_IJK_elem_base_double::set_runge_kutta(int rk_step, double dt_tot, IJK_Field_vector3_double& current_fluxes, IJK_Field_vector3_double& RK3_F_fluxes)
283{
284 if (rk_step == -1)
285 {
287 rk_step_ = -1;
288 dt_tot_ = 0;
289 current_fluxes_ = nullptr;
290 RK3_F_fluxes_ = nullptr;
291 }
292 else
293 {
295 rk_step_ = rk_step;
296 dt_tot_ = dt_tot;
297 current_fluxes_ = &current_fluxes;
298 RK3_F_fluxes_ = &RK3_F_fluxes;
299 }
300}
301
303{
304 compute_(dx, 0);
305}
306
308{
309 compute_(dx, 1);
310}
311
312void Operateur_IJK_elem_base_double::compute_(IJK_Field_double& dx, bool add)
313{
314 IJK_Field_local_double storage;
315 storage.allocate(dx.ni()+1, dx.nj()+1, 4, 0);
316 IJK_Field_local_double tmp[4];
317 for (int i = 0; i < 4; i++)
318 tmp[i].ref_ij(storage, i);
319 IJK_Field_local_double *flux_zmin = &tmp[0];
320 IJK_Field_local_double *flux_zmax = &tmp[1]; // flux a travers la face zmax des volumes de controles
321 IJK_Field_local_double *const flux_x = &tmp[2]; // pointer will not change, content will change
322 IJK_Field_local_double *const flux_y = &tmp[3];
323
324 compute_flux_z(*flux_zmin, 0);
325 correct_flux(flux_zmin, 0, 2);
326 /*
327 * TODO: correct_flux_z_min();
328 * fluxes_to_correct ?
329 * Critere fabs(fluxes_to_correct)<DMAX_FLOAT ?
330 * Zero valeur admissible
331 */
332 const int kmax = dx.nk();
333 for (int k = 0; k < kmax; k++)
334 {
335 compute_flux_x(*flux_x, k);
336 compute_flux_y(*flux_y, k);
337 compute_flux_z(*flux_zmax, k+1);
338 /*
339 * Correct the fluxes before it is summed
340 * with Operator_IJK_div
341 */
342 correct_flux(flux_x, k, 0);
343 correct_flux(flux_y, k, 1);
344 correct_flux(flux_zmax, k+1, 2);
345 Operator_IJK_div(*flux_x, *flux_y, *flux_zmin, *flux_zmax, dx, k, add);
346 swap_data(flux_zmin, flux_zmax); // conserve le flux en z pour la couche z suivante
347 }
348}
349
350void Operateur_IJK_elem_base_double::Operator_IJK_div(const IJK_Field_local_double& flux_x, const IJK_Field_local_double& flux_y,
351 const IJK_Field_local_double& flux_zmin, const IJK_Field_local_double& flux_zmax,
352 IJK_Field_double& resu, int k_layer, bool add)
353{
354 ConstIJK_double_ptr fx(flux_x, 0, 0, 0);
355 ConstIJK_double_ptr fy(flux_y, 0, 0, 0);
356 ConstIJK_double_ptr fzmin(flux_zmin, 0, 0, 0);
357 ConstIJK_double_ptr fzmax(flux_zmax, 0, 0, 0);
358 IJK_double_ptr resu_ptr(resu, 0, 0, k_layer);
359
360
361 const int imax = resu.ni();
362 const int jmax = resu.nj();
363 const int vsize = Simd_double::size();
364 for (int j = 0; ; j++)
365 {
366 for (int i = 0; i < imax; i += vsize)
367 {
368 Simd_double r, a, b;
369 fx.get_center_right(DIRECTION::X, i, a, b);
370 // correct_flux_spherical(a, b, i, j, k_layer, 0);
371 r = a-b;
372
373 fy.get_center_right(DIRECTION::Y,i, a, b);
374 // correct_flux_spherical(a, b, i, j, k_layer, 1);
375 r += a-b;
376
377 fzmin.get_center(i, a);
378 fzmax.get_center(i, b);
379 // correct_flux_spherical(a, b, i, j, k_layer, 2);
380 r += a-b;
381
382 if (add)
383 {
384 resu_ptr.get_center(i, a);
385 r += a;
386 }
387 resu_ptr.put_val(i, r);
388 }
389 // do not execute end_iloop at last iteration (because of assert on valid j+1)
390 if (j+1==jmax)
391 break;
392 fx.next_j();
393 fy.next_j();
394 fzmin.next_j();
395 fzmax.next_j();
396 resu_ptr.next_j();
397 }
398
399}
400
401void Operateur_IJK_elem_base_double::correct_flux(IJK_Field_local_double *const flux, const int k_layer, const int dir)
402{
404 {
405 {
406 int ni = (dir == 0) ? flux->ni() : flux->ni() - 1;
407 int nj = (dir == 1) ? flux->nj() : flux->nj() - 1;
408 for (int j = 0; j < nj; j++)
409 {
410 for (int i = 0; i < ni; i++)
411 {
412 (*current_fluxes_)[dir](i,j,k_layer) = (*flux)(i,j,0);
413 }
414 }
415 }
416
417 runge_kutta3_update_surfacic_fluxes((*current_fluxes_)[dir], (*RK3_F_fluxes_)[dir], rk_step_, k_layer, dt_tot_);
418
419 {
420 int ni = (dir == 0) ? flux->ni() : flux->ni() - 1;
421 int nj = (dir == 1) ? flux->nj() : flux->nj() - 1;
422 for (int j = 0; j < nj; j++)
423 {
424 for (int i = 0; i < ni; i++)
425 {
426 (*flux)(i,j,0) = (*current_fluxes_)[dir](i,j,k_layer);
427 }
428 }
429 }
430 }
431}
432
433void Operateur_IJK_elem_base_double::compute_grad(IJK_Field_vector3_double& dx)
434{
435 IJK_Field_local_double storage;
436 storage.allocate(dx[0].ni()+1, dx[0].nj()+1, 4, 0);
437 IJK_Field_local_double tmp[4];
438 for (int i = 0; i < 4; i++)
439 tmp[i].ref_ij(storage, i);
440 IJK_Field_local_double *flux_zmin = &tmp[0];
441 IJK_Field_local_double *flux_zmax = &tmp[1]; // flux a travers la face zmax des volumes de controles
442 IJK_Field_local_double *const flux_x = &tmp[2]; // pointer will not change, content will change
443 IJK_Field_local_double *const flux_y = &tmp[3];
444
445 compute_flux_z(*flux_zmin, 0);
446
447 const int kmax = dx[0].nk();
448 for (int k = 0; k < kmax; k++)
449 {
450 compute_flux_x(*flux_x, k);
451 compute_flux_y(*flux_y, k);
452 compute_flux_z(*flux_zmax, k+1);
453 fill_grad_field_x_(*flux_x, dx, k);
454 fill_grad_field_y_(*flux_y, dx, k);
455 fill_grad_field_z_(*flux_zmin, *flux_zmax, dx, k);
456 swap_data(flux_zmin, flux_zmax); // conserve le flux en z pour la couche z suivante
457 }
458}
459
461{
462 IJK_Field_local_double storage;
463 storage.allocate(dx.ni()+1, dx.nj()+1, 1, 0);
464 IJK_Field_local_double tmp;
465 tmp.ref_ij(storage, 0);
466 IJK_Field_local_double *const flux_x = &tmp;
467 const int kmax = dx.nk();
468 for (int k = 0; k < kmax; k++)
469 {
470 compute_flux_x(*flux_x, k);
471 fill_grad_field_x_y_(*flux_x, dx, k, 0);
472 }
473}
474
476{
477 IJK_Field_local_double storage;
478 storage.allocate(dx.ni()+1, dx.nj()+1, 1, 0);
479 IJK_Field_local_double tmp;
480 tmp.ref_ij(storage, 0);
481 IJK_Field_local_double *const flux_y = &tmp;
482 const int kmax = dx.nk();
483 for (int k = 0; k < kmax; k++)
484 {
485 compute_flux_y(*flux_y, k);
486 fill_grad_field_x_y_(*flux_y, dx, k, 1);
487 }
488}
489
491{
492 IJK_Field_local_double storage;
493 storage.allocate(dx.ni()+1, dx.nj()+1, 2, 0);
494 IJK_Field_local_double tmp[2];
495 for (int i = 0; i < 2; i++)
496 tmp[i].ref_ij(storage, i);
497 IJK_Field_local_double *flux_zmin = &tmp[0];
498 IJK_Field_local_double *flux_zmax = &tmp[1];
499 compute_flux_z(*flux_zmin, 0);
500 const int kmax = dx.nk();
501 for (int k = 0; k < kmax; k++)
502 {
503 compute_flux_z(*flux_zmax, k+1);
504 fill_grad_field_z_(*flux_zmin, *flux_zmax, dx, k);
505 swap_data(flux_zmin, flux_zmax); // conserve le flux en z pour la couche z suivante
506 }
507}
508
509void Operateur_IJK_elem_base_double::fill_grad_field_x_(IJK_Field_local_double& flux, IJK_Field_vector3_double& resu, int k)
510{
511 fill_grad_field_x_y_(flux, resu[0], k, 0);
512}
513
514void Operateur_IJK_elem_base_double::fill_grad_field_y_(IJK_Field_local_double& flux, IJK_Field_vector3_double& resu, int k)
515{
516 fill_grad_field_x_y_(flux, resu[1], k, 1);
517}
518
519void Operateur_IJK_elem_base_double::fill_grad_field_z_(IJK_Field_local_double& flux_min, IJK_Field_local_double& flux_max, IJK_Field_vector3_double& resu, int k)
520{
521 fill_grad_field_z_(flux_min, flux_max, resu[2], k);
522}
void get_center(int i_offset, _TYPE_ &center) const
returns field(i+i_offset, j, k)
void next_j()
increments the pointer by j_stride (eg, j = j+1)
void get_center_right(DIRECTION _DIR_, int i_offset, _TYPE_ &center, _TYPE_ &right) const
This class encapsulates all the information related to the eulerian mesh for TrioIJK.
Definition Domaine_IJK.h:47
int get_offset_local(int direction) const
Returns the local offset in requested direction.
int get_nb_elem_local(int direction) const
Returns the number of elements owned by this processor in the given direction.
double get_constant_delta(int direction) const
Returns the size of cells in a direction.
const ArrOfDouble & get_delta(int direction) const
Returns the array of mesh cell sizes in requested direction.
Class defining operators and methods for all reading operation in an input flow (file,...
Definition Entree.h:42
void allocate(int ni, int nj, int nk, int ghosts, int additional_k_layers=0, int nb_compo=1, bool external_storage=false)
void ref_ij(IJK_Field_local_template &src, int k_layer)
const Domaine_IJK & get_domaine() const
void put_val(int i_offset, const _TYPE_ &val)
Performs the assignment: field(i+i_offset,j,k) = val.
Definition IJK_ptr.h:32
classe Objet_U Cette classe est la classe de base des Objets de TRUST
Definition Objet_U.h:73
virtual Entree & readOn(Entree &)
Lecture d'un Objet_U sur un flot d'entree Methode a surcharger.
Definition Objet_U.cpp:293
virtual Sortie & printOn(Sortie &) const
Ecriture de l'objet sur un flot de sortie Methode a surcharger.
Definition Objet_U.cpp:282
virtual void Operator_IJK_div(const IJK_Field_local_double &flux_x, const IJK_Field_local_double &flux_y, const IJK_Field_local_double &flux_zmin, const IJK_Field_local_double &flux_zmax, IJK_Field_double &resu, int k_layer, bool add)
virtual void compute_add(IJK_Field_double &dx)
virtual void compute_flux_z(IJK_Field_local_double &resu, const int k_layer)=0
virtual void compute_grad(IJK_Field_vector3_double &dx)
virtual void compute_set(IJK_Field_double &dx)
virtual void correct_flux(IJK_Field_local_double *const flux, const int k_layer, const int dir)
virtual void compute_flux_y(IJK_Field_local_double &resu, const int k_layer)=0
virtual void compute_flux_x(IJK_Field_local_double &resu, const int k_layer)=0
virtual void compute_grad_x(IJK_Field_double &dx)
virtual void compute_grad_y(IJK_Field_double &dx)
void set_runge_kutta(int rk_step, double dt_tot, IJK_Field_vector3_double &current_fluxes, IJK_Field_vector3_double &RK3_F_fluxes)
IJK_Field_vector3_double * current_fluxes_
virtual void compute_grad_z(IJK_Field_double &dx)
IJK_Field_vector3_double * RK3_F_fluxes_
virtual void compute_flux_y_vy(IJK_Field_local_double &resu, const int k_layer)=0
virtual void compute_flux_y_vx(IJK_Field_local_double &resu, const int k_layer)=0
virtual void exec_after_divergence_flux_y(IJK_Field_double &resu, const int k_layer)
virtual void compute_flux_x_vx(IJK_Field_local_double &resu, const int k_layer)=0
virtual void exec_after_divergence_flux_x(IJK_Field_double &resu, const int k_layer)
virtual void compute_flux_z_vz(IJK_Field_local_double &resu, const int k_layer)=0
virtual void correct_flux(IJK_Field_local_double *const flux, const int k_layer, const int dir, const int compo)
void compute_set(IJK_Field_double &dvx, IJK_Field_double &dvy, IJK_Field_double &dvz)
virtual double compute_dtstab_convection_local(IJK_Field_double &dvx, IJK_Field_double &dvy, IJK_Field_double &dvz)
virtual void compute_flux_x_vz(IJK_Field_local_double &resu, const int k_layer)=0
virtual void compute_flux_y_vz(IJK_Field_local_double &resu, const int k_layer)=0
virtual void compute_flux_x_vy(IJK_Field_local_double &resu, const int k_layer)=0
virtual void compute_flux_z_vy(IJK_Field_local_double &resu, const int k_layer)=0
virtual void compute_flux_z_vx(IJK_Field_local_double &resu, const int k_layer)=0
void compute_add(IJK_Field_double &dvx, IJK_Field_double &dvy, IJK_Field_double &dvz)
virtual void exec_after_divergence_flux_z(IJK_Field_double &resu, const int k_layer)
virtual void Operator_IJK_div(const IJK_Field_local_double &flux_x, const IJK_Field_local_double &flux_y, const IJK_Field_local_double &flux_zmin, const IJK_Field_local_double &flux_zmax, IJK_Field_double &resu, int k_layer, bool add)
Classe de base des flux de sortie.
Definition Sortie.h:52