Flow123d  release_2.1.0-87-gfbc1563
transport.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file transport.cc
15  * @ingroup transport
16  * @brief Transport
17  */
18 
19 #include <memory>
20 
21 #include "system/system.hh"
22 #include "system/sys_profiler.hh"
23 
24 #include "mesh/mesh.h"
25 #include "mesh/partitioning.hh"
26 #include "transport/transport.h"
27 
28 #include "la/distribution.hh"
29 
30 #include "la/sparse_graph.hh"
31 #include <iostream>
32 #include <iomanip>
33 #include <string>
34 
35 #include "io/output_time.hh"
36 #include "tools/time_governor.hh"
37 #include "coupling/balance.hh"
38 #include "input/accessors.hh"
39 #include "input/input_type.hh"
40 
42 #include "fields/field_values.hh"
44 #include "fields/generic_field.hh"
45 
46 #include "reaction/isotherm.hh" // SorptionType enum
47 
48 
49 FLOW123D_FORCE_LINK_IN_CHILD(convectionTransport);
50 
51 
52 namespace IT = Input::Type;
53 
54 const string _equation_name = "Solute_Advection_FV";
55 
57  Input::register_class< ConvectionTransport, Mesh &, const Input::Record >(_equation_name) +
59 
61 {
62  return IT::Record(_equation_name, "Explicit in time finite volume method for advection only solute transport.")
64  .declare_key("input_fields", IT::Array(
65  EqData().make_field_descriptor_type(_equation_name)),
67  "")
68  .declare_key("output",
69  EqData().output_fields.make_output_type(_equation_name, ""),
70  IT::Default("{ \"fields\": [ \"conc\" ] }"),
71  "Setting of the fields output.")
72  .close();
73 }
74 
75 
77 {
78  ADD_FIELD(bc_conc, "Boundary conditions for concentrations.", "0.0");
79  bc_conc.units( UnitSI().kg().m(-3) );
80  ADD_FIELD(init_conc, "Initial concentrations.", "0.0");
81  init_conc.units( UnitSI().kg().m(-3) );
82 
83  output_fields += *this;
84  output_fields += conc_mobile.name("conc")
85  .units( UnitSI().kg().m(-3) )
87  output_fields += region_id.name("region_id")
90 }
91 
92 
94 : ConcentrationTransportBase(init_mesh, in_rec),
95  is_mass_diag_changed(false),
96  input_rec(in_rec),
97  mh_dh(nullptr),
98  sources_corr(nullptr)
99 {
100  START_TIMER("ConvectionTransport");
101  this->eq_data_ = &data_;
102 
103  transport_matrix_time = -1.0; // or -infty
104  transport_bc_time = -1.0;
106  is_src_term_scaled = false;
107  is_bc_term_scaled = false;
108 }
109 
111 {
113 
115 
117  data_.set_input_list( input_rec.val<Input::Array>("input_fields") );
118  data_.set_mesh(*mesh_);
119 
123 
124  // register output vectors
130  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
131  {
132  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
133  auto output_field_ptr = out_conc[sbi].create_field<3, FieldValue<3>::Scalar>(n_substances());
134  data_.conc_mobile[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
135  }
136  //output_stream_->add_admissible_field_names(input_rec.val<Input::Array>("output_fields"));
137  //output_stream_->mark_output_times(*time_);
138 
140  //cout << "Transport." << endl;
141  //cout << time().marks();
142 
143  if (balance_ != nullptr)
144  balance_->allocate(el_ds->lsize(), 1);
145 
146 }
147 
148 
149 //=============================================================================
150 // MAKE TRANSPORT
151 //=============================================================================
153 
154 // int * id_4_old = new int[mesh_->n_elements()];
155 // int i = 0;
156 // FOR_ELEMENTS(mesh_, ele) id_4_old[i++] = ele.index();
157 // mesh_->get_part()->id_maps(mesh_->n_elements(), id_4_old, el_ds, el_4_loc, row_4_el);
158 // delete[] id_4_old;
159  el_ds = mesh_->get_el_ds();
162 
163  // TODO: make output of partitioning is usefull but makes outputs different
164  // on different number of processors, which breaks tests.
165  //
166  // Possible solution:
167  // - have flag in ini file to turn this output ON
168  // - possibility to have different ref_output for different num of proc.
169  // - or do not test such kind of output
170  //
171  //FOR_ELEMENTS(mesh_, ele) {
172  // ele->pid=el_ds->get_proc(row_4_el[ele.index()]);
173  //}
174 
175 }
176 
177 
178 
180 {
181  unsigned int sbi;
182 
183  if (sources_corr) {
184  //Destroy mpi vectors at first
185  MatDestroy(&tm);
186  VecDestroy(&mass_diag);
187  VecDestroy(&vpmass_diag);
188  VecDestroy(&vcfl_flow_);
189  VecDestroy(&vcfl_source_);
190  delete cfl_flow_;
191  delete cfl_source_;
192 
193  for (sbi = 0; sbi < n_substances(); sbi++) {
194  // mpi vectors
195  VecDestroy(&(vconc[sbi]));
196  VecDestroy(&(vpconc[sbi]));
197  VecDestroy(&(bcvcorr[sbi]));
198  VecDestroy(&(vcumulative_corr[sbi]));
199  VecDestroy(&(v_tm_diag[sbi]));
200  VecDestroy(&(v_sources_corr[sbi]));
201 
202  // arrays of arrays
203  delete conc[sbi];
204  delete cumulative_corr[sbi];
205  delete tm_diag[sbi];
206  delete sources_corr[sbi];
207  }
208 
209  // arrays of mpi vectors
210  delete vconc;
211  delete vpconc;
212  delete bcvcorr;
213  delete vcumulative_corr;
214  delete v_tm_diag;
215  delete v_sources_corr;
216 
217  // arrays of arrays
218  delete conc;
219  delete cumulative_corr;
220  delete tm_diag;
221  delete sources_corr;
222  }
223 }
224 
225 
226 
227 
228 
230 {
231  FOR_ELEMENTS(mesh_, elem)
232  {
233  if (!el_ds->is_local(row_4_el[elem.index()])) continue;
234 
235  unsigned int index = row_4_el[elem.index()] - el_ds->begin();
236  ElementAccessor<3> ele_acc = mesh_->element_accessor(elem.index());
237  arma::vec value = data_.init_conc.value(elem->centre(), ele_acc);
238 
239  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
240  conc[sbi][index] = value(sbi);
241  }
242 
243 }
244 
245 //=============================================================================
246 // ALLOCATE OF TRANSPORT VARIABLES (ELEMENT & NODES)
247 //=============================================================================
249 
250  unsigned int i, sbi, n_subst;
251  n_subst = n_substances();
252 
253  sources_corr = new double*[n_subst];
254  tm_diag = new double*[n_subst];
255  cumulative_corr = new double*[n_subst];
256  for (sbi = 0; sbi < n_subst; sbi++) {
257  cumulative_corr[sbi] = new double[el_ds->lsize()];
258  sources_corr[sbi] = new double[el_ds->lsize()];
259  tm_diag[sbi] = new double[el_ds->lsize()];
260  }
261 
262  conc = new double*[n_subst];
263  out_conc.clear();
264  out_conc.resize(n_subst);
265  for (sbi = 0; sbi < n_subst; sbi++) {
266  conc[sbi] = new double[el_ds->lsize()];
267  out_conc[sbi].resize( el_ds->size() );
268  for (i = 0; i < el_ds->lsize(); i++) {
269  conc[sbi][i] = 0.0;
270  }
271  }
272 
273  cfl_flow_ = new double[el_ds->lsize()];
274  cfl_source_ = new double[el_ds->lsize()];
275 }
276 
277 //=============================================================================
278 // ALLOCATION OF TRANSPORT VECTORS (MPI)
279 //=============================================================================
281 
282  int sbi, n_subst, rank, np;
283  n_subst = n_substances();
284 
285  MPI_Barrier(PETSC_COMM_WORLD);
286  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
287  MPI_Comm_size(PETSC_COMM_WORLD, &np);
288 
289  vconc = new Vec[n_subst];
290  vpconc = new Vec[n_subst];
291  bcvcorr = new Vec[n_subst];
292  vcumulative_corr = new Vec[n_subst];
293  v_tm_diag = new Vec[n_subst];
294  v_sources_corr = new Vec[n_subst];
295 
296 
297  for (sbi = 0; sbi < n_subst; sbi++) {
298  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &bcvcorr[sbi]);
299  VecZeroEntries(bcvcorr[sbi]);
300  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(), conc[sbi],
301  &vconc[sbi]);
302 
303  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpconc[sbi]);
304  VecZeroEntries(vconc[sbi]);
305  VecZeroEntries(vpconc[sbi]);
306 
307  // SOURCES
308  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
309  cumulative_corr[sbi],&vcumulative_corr[sbi]);
310 
311  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
312  sources_corr[sbi],&v_sources_corr[sbi]);
313 
314  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
315  tm_diag[sbi],&v_tm_diag[sbi]);
316 
317  VecZeroEntries(vcumulative_corr[sbi]);
318  VecZeroEntries(out_conc[sbi].get_data_petsc());
319  }
320 
321 
322  MatCreateAIJ(PETSC_COMM_WORLD, el_ds->lsize(), el_ds->lsize(), mesh_->n_elements(),
323  mesh_->n_elements(), 16, PETSC_NULL, 4, PETSC_NULL, &tm);
324 
325  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &mass_diag);
326  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpmass_diag);
327 
328  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
330  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
332 }
333 
334 
336 {
337  START_TIMER ("set_boundary_conditions");
338 
340 
341  unsigned int sbi, loc_el, loc_b = 0;
342 
343  // Assembly bcvcorr vector
344  for(sbi=0; sbi < n_substances(); sbi++) VecZeroEntries(bcvcorr[sbi]);
345 
346  if (balance_ != nullptr)
347  balance_->start_flux_assembly(subst_idx);
348 
349  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
350  elm = mesh_->element(el_4_loc[loc_el]);
351  if (elm->boundary_idx_ != NULL) {
352  unsigned int new_i = row_4_el[elm.index()];
353 
354  FOR_ELEMENT_SIDES(elm,si) {
355  Boundary *b = elm->side(si)->cond();
356  if (b != NULL) {
357  double flux = mh_dh->side_flux( *(elm->side(si)) );
358  if (flux < 0.0) {
359  double aij = -(flux / elm->measure() );
360 
361  arma::vec value = data_.bc_conc.value( b->element()->centre(), b->element_accessor() );
362  for (sbi=0; sbi<n_substances(); sbi++)
363  VecSetValue(bcvcorr[sbi], new_i, value[sbi] * aij, ADD_VALUES);
364 
365  if (balance_ != nullptr)
366  {
367  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
368  {
369  // CAUTION: It seems that PETSc possibly optimize allocated space during assembly.
370  // So we have to add also values that may be non-zero in future due to changing velocity field.
371  balance_->add_flux_matrix_values(subst_idx[sbi], loc_b, {row_4_el[el_4_loc[loc_el]]}, {0.});
372  balance_->add_flux_vec_value(subst_idx[sbi], loc_b, flux*value[sbi]);
373  }
374  }
375  } else {
376  for (sbi=0; sbi<n_substances(); sbi++)
377  VecSetValue(bcvcorr[sbi], new_i, 0, ADD_VALUES);
378 
379  if (balance_ != nullptr)
380  {
381  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
382  {
383  balance_->add_flux_matrix_values(subst_idx[sbi], loc_b, {row_4_el[el_4_loc[loc_el]]}, {flux});
384  balance_->add_flux_vec_value(subst_idx[sbi], loc_b, 0);
385  }
386  }
387  }
388  ++loc_b;
389  }
390  }
391 
392  }
393  }
394 
395  if (balance_ != nullptr)
396  balance_->finish_flux_assembly(subst_idx);
397 
398  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyBegin(bcvcorr[sbi]);
399  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyEnd(bcvcorr[sbi]);
400 
401  // we are calling set_boundary_conditions() after next_time() and
402  // we are using data from t() before, so we need to set corresponding bc time
404 }
405 
406 
407 //=============================================================================
408 // COMPUTE SOURCES
409 //=============================================================================
411 
412  //temporary variables
413  unsigned int loc_el, sbi;
414  double csection, source, diag;
415 
416  Element *ele;
417  ElementAccessor<3> ele_acc;
418  arma::vec3 p;
419  arma::vec src_density(n_substances()), src_conc(n_substances()), src_sigma(n_substances());
420 
421  //TODO: would it be possible to check the change in data for chosen substance? (may be in multifields?)
422 
423  //checking if the data were changed
424  if( (data_.sources_density.changed() )
425  || (data_.sources_conc.changed() )
426  || (data_.sources_sigma.changed() )
427  || (data_.cross_section.changed()))
428  {
429  START_TIMER("sources_reinit");
430  if (balance_ != nullptr) balance_->start_source_assembly(subst_idx);
431 
432  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++)
433  {
434  ele = mesh_->element(el_4_loc[loc_el]);
435  ele_acc = ele->element_accessor();
436  p = ele_acc.centre();
437  csection = data_.cross_section.value(p, ele_acc);
438 
439  // read for all substances
440  src_density = data_.sources_density.value(p, ele_acc);
441  src_conc = data_.sources_conc.value(p, ele_acc);
442  src_sigma = data_.sources_sigma.value(p, ele_acc);
443 
444  double max_cfl=0;
445  for (sbi = 0; sbi < n_substances(); sbi++)
446  {
447  source = csection * (src_density(sbi) + src_sigma(sbi) * src_conc(sbi));
448  // addition to RHS
449  sources_corr[sbi][loc_el] = source;
450  // addition to diagonal of the transport matrix
451  diag = src_sigma(sbi) * csection;
452  tm_diag[sbi][loc_el] = - diag;
453 
454  // compute maximal cfl condition over all substances
455  max_cfl = std::max(max_cfl, fabs(diag));
456 
457  if (balance_ != nullptr)
458  {
459  balance_->add_source_matrix_values(sbi, ele_acc.region().bulk_idx(), {row_4_el[el_4_loc[loc_el]]},
460  {- src_sigma(sbi) * ele->measure() * csection});
461  balance_->add_source_vec_values(sbi, ele_acc.region().bulk_idx(), {row_4_el[el_4_loc[loc_el]]},
462  {source * ele->measure()});
463  }
464  }
465 
466  cfl_source_[loc_el] = max_cfl;
467  }
468 
469  if (balance_ != nullptr) balance_->finish_source_assembly(subst_idx);
470 
471  END_TIMER("sources_reinit");
472  }
473 }
474 
475 
476 
478 {
480 
483 
486 
487  if (balance_ != nullptr)
488  {
489  START_TIMER("Convection balance zero time step");
490 
494 
496  }
497 
498  // write initial condition
499  output_data();
500 }
501 
502 
504 {
505  OLD_ASSERT(mh_dh, "Null MH object.\n" );
506  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
507 
508  START_TIMER("data reinit");
509 
510  bool cfl_changed = false;
511 
512  // if FLOW or DATA changed ---------------------> recompute transport matrix
514  {
517  cfl_changed = true;
518  DebugOut() << "CFL changed - flow.\n";
519  }
520 
522  {
524  cfl_changed = true;
525  DebugOut() << "CFL changed - mass matrix.\n";
526  }
527 
528  // if DATA changed ---------------------> recompute concentration sources (rhs and matrix diagonal)
531  {
533  is_src_term_scaled = false;
534  cfl_changed = true;
535  DebugOut() << "CFL changed - source.\n";
536  }
537 
538  // now resolve the CFL condition
539  if(cfl_changed)
540  {
541  // find maximum of sum of contribution from flow and sources: MAX(vcfl_flow_ + vcfl_source_)
542  Vec cfl;
543  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(),PETSC_DETERMINE, &cfl);
544  VecWAXPY(cfl, 1.0, vcfl_flow_, vcfl_source_);
545  VecMaxPointwiseDivide(cfl,mass_diag, &cfl_max_step);
546  // get a reciprocal value as a time constraint
548  DebugOut().fmt("CFL constraint (transport): {}\n", cfl_max_step);
549  }
550 
551  // although it does not influence CFL, compute BC so the full system is assembled
553  || data_.porosity.changed()
555  || data_.bc_conc.changed() )
556  {
558  is_bc_term_scaled = false;
559  }
560 
561  END_TIMER("data reinit");
562 
563  // return time constraint
564  time_constraint = cfl_max_step;
565  return cfl_changed;
566 }
567 
569 
570  START_TIMER("convection-one step");
571 
572  // proceed to next time - which we are about to compute
573  // explicit scheme looks one step back and uses data from previous time
574  // (data time set previously in assess_time_constraint())
575  time_->next_time();
576 
577  double dt_new = time_->dt(), // current time step we are about to compute
578  dt_scaled = dt_new / time_->last_dt(); // scaling ratio to previous time step
579 
580  START_TIMER("time step rescaling");
581 
582  // if FLOW or DATA or BC or DT changed ---------------------> rescale boundary condition
584  {
585  DebugOut() << "BC - rescale dt.\n";
586  //choose between fresh scaling with new dt or rescaling to a new dt
587  double dt = (!is_bc_term_scaled) ? dt_new : dt_scaled;
588  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
589  VecScale(bcvcorr[sbi], dt);
590  is_bc_term_scaled = true;
591  }
592 
593 
594  // if DATA or TIME STEP changed -----------------------> rescale source term
596  DebugOut() << "SRC - rescale dt.\n";
597  //choose between fresh scaling with new dt or rescaling to a new dt
598  double dt = (!is_src_term_scaled) ? dt_new : dt_scaled;
599  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
600  {
601  VecScale(v_sources_corr[sbi], dt);
602  VecScale(v_tm_diag[sbi], dt);
603  }
604  is_src_term_scaled = true;
605  }
606 
607  // if DATA or TIME STEP changed -----------------------> rescale transport matrix
609  DebugOut() << "TM - rescale dt.\n";
610  //choose between fresh scaling with new dt or rescaling to a new dt
611  double dt = (!is_convection_matrix_scaled) ? dt_new : dt_scaled;
612 
613  MatScale(tm, dt);
615  }
616 
617  END_TIMER("time step rescaling");
618 
619 
620  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
622  {
623  VecCopy(mass_diag, vpmass_diag);
625  } else is_mass_diag_changed = false;
626 
627 
628  // Compute new concentrations for every substance.
629 
630  for (unsigned int sbi = 0; sbi < n_substances(); sbi++) {
631  // one step in MOBILE phase
632  START_TIMER("mat mult");
633 
634  // tm_diag is a diagonal part of transport matrix, which depends on substance data (sources_sigma)
635  // Wwe need keep transport matrix independent of substance, therefore we keep this diagonal part
636  // separately in a vector tm_diag.
637  // Computation: first, we compute this diagonal addition D*pconc and save it temporaly into RHS
638 
639  // RHS = D*pconc, where D is diagonal matrix represented by a vector
640  VecPointwiseMult(vcumulative_corr[sbi], v_tm_diag[sbi], vconc[sbi]); //w = x.*y
641 
642  // Then we add boundary terms ans other source terms into RHS.
643  // RHS = 1.0 * bcvcorr + 1.0 * v_sources_corr + 1.0 * rhs
644  VecAXPBYPCZ(vcumulative_corr[sbi], 1.0, 1.0, 1.0, bcvcorr[sbi], v_sources_corr[sbi]); //z = ax + by + cz
645 
646  // Then we set the new previous concentration.
647  VecCopy(vconc[sbi], vpconc[sbi]); // pconc = conc
648  // And finally proceed with transport matrix multiplication.
649  if (is_mass_diag_changed) {
650  VecPointwiseMult(vconc[sbi], vconc[sbi], vpmass_diag); // vconc*=vpmass_diag
651  MatMultAdd(tm, vpconc[sbi], vconc[sbi], vconc[sbi]); // vconc+=tm*vpconc
652  VecAXPY(vconc[sbi], 1, vcumulative_corr[sbi]); // vconc+=vcumulative_corr
653  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
654  } else {
655  MatMultAdd(tm, vpconc[sbi], vcumulative_corr[sbi], vconc[sbi]); // vconc =tm*vpconc+vcumulative_corr
656  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
657  VecAXPY(vconc[sbi], 1, vpconc[sbi]); // vconc+=vpconc
658  }
659 
660  END_TIMER("mat mult");
661  }
662 
663  if (balance_ != nullptr)
665 
666  END_TIMER("convection-one step");
667 }
668 
669 
670 void ConvectionTransport::set_target_time(double target_time)
671 {
672 
673  //sets target_mark_type (it is fixed) to be met in next_time()
674  time_->marks().add(TimeMark(target_time, target_mark_type));
675 
676  // This is done every time TOS calls update_solution.
677  // If CFL condition is changed, time fixation will change later from TOS.
678 
679  // Set the same constraint as was set last time.
680 
681  // TODO: fix this hack, remove this method completely, leaving just the first line at the calling point
682  // in TransportOperatorSplitting::update_solution()
683  // doing this directly leads to choose of large time step violationg CFL condition
684  if (cfl_max_step > time_->dt()*1e-10)
685  time_->set_upper_constraint(cfl_max_step, "CFL condition used from previous step.");
686 
687  // fixing convection time governor till next target_mark_type (got from TOS or other)
688  // may have marks for data changes
690 }
691 
692 
694 {
696 
697  VecZeroEntries(mass_diag);
698 
699  if (balance_ != nullptr)
700  balance_->start_mass_assembly(subst_idx);
701 
702  for (unsigned int loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
703  elm = mesh_->element(el_4_loc[loc_el]);
704 
705  double csection = data_.cross_section.value(elm->centre(), elm->element_accessor());
706  //double por_m = data_.porosity.value(elm->centre(), elm->element_accessor());
707  double por_m = data_.water_content.value(elm->centre(), elm->element_accessor());
708 
709  if (balance_ != nullptr)
710  {
711  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
712  balance_->add_mass_matrix_values(subst_idx[sbi], elm->region().bulk_idx(), {row_4_el[el_4_loc[loc_el]]}, {csection*por_m*elm->measure()} );
713  }
714 
715  VecSetValue(mass_diag, row_4_el[elm.index()], csection*por_m, INSERT_VALUES);
716  }
717 
718  if (balance_ != nullptr)
719  balance_->finish_mass_assembly(subst_idx);
720 
721  VecAssemblyBegin(mass_diag);
722  VecAssemblyEnd(mass_diag);
723 
724  is_mass_diag_changed = true;
725 }
726 
727 
728 //=============================================================================
729 // CREATE TRANSPORT MATRIX
730 //=============================================================================
732 
733  START_TIMER("convection_matrix_assembly");
734 
737  struct Edge *edg;
738  unsigned int n;
739  int s, j, np, rank, new_j, new_i;
740  double aij, aii;
741 
742  MatZeroEntries(tm);
743 
744  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
745  MPI_Comm_size(PETSC_COMM_WORLD, &np);
746 
747  double flux, flux2, edg_flux;
748 
749  aii = 0.0;
750 
751  for (unsigned int loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
752  elm = mesh_->element(el_4_loc[loc_el]);
753  new_i = row_4_el[elm.index()];
754 
755  FOR_ELEMENT_SIDES(elm,si) {
756  // same dim
757  flux = mh_dh->side_flux( *(elm->side(si)) );
758  if (elm->side(si)->cond() == NULL) {
759  edg = elm->side(si)->edge();
760  edg_flux = 0;
761  for( int s=0; s < edg->n_sides; s++) {
762  flux2 = mh_dh->side_flux( *(edg->side(s)) );
763  if ( flux2 > 0) edg_flux+= flux2;
764  }
765  FOR_EDGE_SIDES(edg,s)
766  // this test should also eliminate sides facing to lower dim. elements in comp. neighboring
767  // These edges on these sides should have just one side
768  if (edg->side(s) != elm->side(si)) {
769  j = ELEMENT_FULL_ITER(mesh_, edg->side(s)->element()).index();
770  new_j = row_4_el[j];
771 
772  flux2 = mh_dh->side_flux( *(edg->side(s)));
773  if ( flux2 > 0.0 && flux <0.0)
774  aij = -(flux * flux2 / ( edg_flux * elm->measure() ) );
775  else aij =0;
776  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
777  }
778  }
779  if (flux > 0.0)
780  aii -= (flux / elm->measure() );
781  } // end same dim //ELEMENT_SIDES
782 
783  FOR_ELM_NEIGHS_VB(elm,n) // comp model
784  {
785  el2 = ELEMENT_FULL_ITER(mesh_, elm->neigh_vb[n]->side()->element() ); // higher dim. el.
786  OLD_ASSERT( el2 != elm, "Elm. same\n");
787  new_j = row_4_el[el2.index()];
788  flux = mh_dh->side_flux( *(elm->neigh_vb[n]->side()) );
789 
790  // volume source - out-flow from higher dimension
791  if (flux > 0.0) aij = flux / elm->measure();
792  else aij=0;
793  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
794  // out flow from higher dim. already accounted
795 
796  // volume drain - in-flow to higher dimension
797  if (flux < 0.0) {
798  aii -= (-flux) / elm->measure(); // diagonal drain
799  aij = (-flux) / el2->measure();
800  } else aij=0;
801  MatSetValue(tm, new_j, new_i, aij, INSERT_VALUES);
802  }
803 
804  MatSetValue(tm, new_i, new_i, aii, INSERT_VALUES);
805 
806  cfl_flow_[loc_el] = fabs(aii);
807  aii = 0.0;
808  } // END ELEMENTS
809 
810  MatAssemblyBegin(tm, MAT_FINAL_ASSEMBLY);
811  MatAssemblyEnd(tm, MAT_FINAL_ASSEMBLY);
812 
814  END_TIMER("convection_matrix_assembly");
815 
817 }
818 
819 
820 
821 
822 
823 //=============================================================================
824 // OUTPUT VECTOR GATHER
825 //=============================================================================
827 
828  unsigned int sbi;
829  IS is;
830 
831  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el, PETSC_COPY_VALUES, &is); //WithArray
832  VecScatterCreate(vconc[0], is, out_conc[0].get_data_petsc(), PETSC_NULL, &vconc_out_scatter);
833  for (sbi = 0; sbi < n_substances(); sbi++) {
834  VecScatterBegin(vconc_out_scatter, vconc[sbi], out_conc[sbi].get_data_petsc(), INSERT_VALUES, SCATTER_FORWARD);
835  VecScatterEnd(vconc_out_scatter, vconc[sbi], out_conc[sbi].get_data_petsc(), INSERT_VALUES, SCATTER_FORWARD);
836  }
837  VecScatterDestroy(&(vconc_out_scatter));
838  ISDestroy(&(is));
839 }
840 
841 
843  return conc;
844 }
845 
846 void ConvectionTransport::get_par_info(int * &el_4_loc_out, Distribution * &el_distribution_out){
847  el_4_loc_out = this->el_4_loc;
848  el_distribution_out = this->el_ds;
849  return;
850 }
851 
852 //int *ConvectionTransport::get_el_4_loc(){
853 // return el_4_loc;
854 //}
855 
857  return row_4_el;
858 }
859 
860 
861 
863 {
864  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
865  {
866  balance_->calculate_cumulative_sources(sbi, vpconc[sbi], time_->dt());
867  balance_->calculate_cumulative_fluxes(sbi, vpconc[sbi], time_->dt());
868  }
869 }
870 
871 
873 {
874  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
875  {
876  balance_->calculate_mass(sbi, vconc[sbi]);
877  balance_->calculate_source(sbi, vconc[sbi]);
878  balance_->calculate_flux(sbi, vconc[sbi]);
879  }
880 }
881 
882 
884 
886  //if ( data_.output_fields.is_field_output_time(data_.conc_mobile, time().step()) ) {
888  //}
889 
890  data_.output_fields.output(time().step());
891 }
892 
893 void ConvectionTransport::set_balance_object(std::shared_ptr<Balance> balance)
894 {
895  balance_ = balance;
896  subst_idx = balance_->add_quantities(substances_.names());
897 }
TimeGovernor & time()
Definition: equation.hh:148
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:201
unsigned int size() const
get global size
void set_input_list(Input::Array input_list)
Definition: field_set.hh:180
void calculate_instant_balance() override
Definition: transport.cc:872
FieldSet * eq_data_
Definition: equation.hh:232
std::vector< VectorSeqDouble > out_conc
Definition: transport.h:305
bool is_mass_diag_changed
Flag indicates that porosity or cross_section changed during last time.
Definition: transport.h:266
double measure() const
Definition: elements.cc:89
double time_changed() const
void update_solution() override
Definition: transport.cc:568
Accessor to input data conforming to declared Array.
Definition: accessors.hh:561
double end_time() const
End time.
double transport_matrix_time
Definition: transport.h:291
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:593
int tlevel() const
const std::vector< std::string > & names()
Definition: substance.hh:85
#define FOR_EDGE_SIDES(i, j)
Definition: edges.h:41
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:56
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:405
void alloc_transport_vectors()
Definition: transport.cc:248
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:200
void output(TimeStep step)
double fix_dt_until_mark()
Fixing time step until fixed time mark.
MultiField< 3, FieldValue< 3 >::Scalar > conc_mobile
Calculated concentrations in the mobile zone.
Definition: transport.h:89
void create_mass_matrix()
Definition: transport.cc:693
void next_time()
Proceed to the next time according to current estimated time step.
void initialize() override
Definition: transport.cc:110
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:310
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:273
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:105
void set_initial_condition()
Definition: transport.cc:229
double transport_bc_time
Time of the last update of the boundary condition terms.
Definition: transport.h:292
double * cfl_flow_
Definition: transport.h:278
MultiField< 3, FieldValue< 3 >::Scalar > sources_density
Concentration sources - density of substance source, only positive part is used.
void set_boundary_conditions()
Definition: transport.cc:335
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
#define ELEMENT_FULL_ITER(_mesh_, i)
Definition: mesh.h:413
Definition: mesh.h:95
Fields computed from the mesh data.
EquationOutput output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
Definition: transport.h:93
void initialize(std::shared_ptr< OutputTime > stream, Input::Record in_rec, const TimeGovernor &tg)
double ** sources_corr
Definition: transport.h:269
#define FOR_ELEMENT_SIDES(i, j)
Definition: elements.h:189
MultiField< 3, FieldValue< 3 >::Scalar > sources_conc
int index() const
Definition: sys_vector.hh:78
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_conc
Definition: transport.h:83
int * get_el_4_loc() const
Definition: mesh.h:162
int n_sides
Definition: edges.h:36
Definition: edges.h:26
const RegionDB & region_db() const
Definition: mesh.h:147
void set_balance_object(std::shared_ptr< Balance > balance) override
Definition: transport.cc:893
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:703
double ** conc
Concentrations for phase, substance, element.
Definition: transport.h:297
const TimeStep & step(int index=-1) const
virtual ~ConvectionTransport()
Definition: transport.cc:179
static Input::Type::Abstract & get_input_type()
Common specification of the input record for secondary equations.
MultiField< 3, FieldValue< 3 >::Scalar > init_conc
Initial concentrations.
Definition: transport.h:86
double t() const
#define ADD_FIELD(name,...)
Definition: field_set.hh:269
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
void zero_time_step() override
Definition: transport.cc:477
static TimeMarks & marks()
void output_vector_gather()
Definition: transport.cc:826
Basic time management class.
Class for declaration of inputs sequences.
Definition: type_base.hh:339
const Input::Record input_rec
Record with input specification.
Definition: transport.h:308
static constexpr bool value
Definition: json.hpp:87
const string _equation_name
Definition: transport.cc:54
double last_t() const
const MH_DofHandler * mh_dh
Definition: transport.h:325
static const Input::Type::Record & get_input_type()
Definition: transport.cc:60
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:193
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void setup_components()
unsigned int n_elements() const
Definition: mesh.h:133
bool is_local(unsigned int idx) const
identify local index
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
TimeMark::Type equation_fixed_mark_type() const
MultiField< 3, FieldValue< 3 >::Scalar > sources_sigma
Concentration sources - Robin type, in_flux = sources_sigma * (sources_conc - mobile_conc) ...
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:58
unsigned int begin(int proc) const
get starting local index
Accessor to the data with type Type::Record.
Definition: accessors.hh:286
const Ret val(const string &key) const
static auto region_id(Mesh &mesh) -> IndexField
ConvectionTransport(Mesh &init_mesh, const Input::Record in_rec)
Definition: transport.cc:93
static const int registrar
Registrar of class to factory.
Definition: transport.h:250
bool is_changed_dt() const
#define START_TIMER(tag)
Starts a timer with specified tag.
SubstanceList substances_
Transported substances.
Definition: transport.h:318
Mesh * mesh_
Definition: equation.hh:223
Element * element()
Definition: boundaries.cc:39
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:350
double ** get_concentration_matrix() override
Getter for array of concentrations per element.
Definition: transport.cc:842
Field< 3, FieldValue< 3 >::Scalar > porosity
Mobile porosity - usually saturated water content in the case of unsaturated flow model...
std::shared_ptr< Balance > balance_
object for calculation and writing the mass balance to file.
Definition: equation.hh:235
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:488
void calculate_cumulative_balance() override
Definition: transport.cc:862
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
Definition: transport.h:328
Distribution * get_el_ds() const
Definition: mesh.h:156
#define MPI_Comm_size
Definition: mpi.h:235
Field< 3, FieldValue< 3 >::Scalar > water_content
Water content - result of unsaturated water flow model or porosity.
Vec * vconc
Concentration vectors for mobile phase.
Definition: transport.h:295
Region region() const
Definition: accessors.hh:98
#define MPI_Comm_rank
Definition: mpi.h:236
void mark_input_times(const TimeGovernor &tg)
Definition: field_set.hh:208
TimeMark add(const TimeMark &mark)
Definition: time_marks.cc:82
Field< 3, FieldValue< 3 >::Integer > region_id
Definition: transport.h:88
void create_transport_matrix_mpi()
Definition: transport.cc:731
Support classes for parallel programing.
void alloc_transport_structs_mpi()
Definition: transport.cc:280
int set_upper_constraint(double upper, std::string message)
Sets upper constraint for the next time step estimating.
double side_flux(const Side &side) const
temporary replacement for DofHandler accessor, flux through given side
bool is_convection_matrix_scaled
Definition: transport.h:263
FLOW123D_FORCE_LINK_IN_CHILD(convectionTransport)
double * cfl_source_
Definition: transport.h:278
#define ELEMENT_FULL_ITER_NULL(_mesh_)
Definition: mesh.h:419
virtual MultiFieldValue::return_type value(const Point &p, const ElementAccessor< spacedim > &elm) const
void set_components(const std::vector< string > &names)
Definition: field_set.hh:167
ElementFullIter element() const
Definition: side_impl.hh:52
void compute_concentration_sources()
Assembles concentration sources for each substance. note: the source of concentration is multiplied b...
Definition: transport.cc:410
bool evaluate_time_constraint(double &time_constraint) override
Definition: transport.cc:503
VecScatter vconc_out_scatter
Definition: transport.h:281
double dt() const
Vec vcfl_flow_
Parallel vector for flow contribution to CFL condition.
Definition: transport.h:276
arma::vec3 centre() const
Definition: elements.cc:120
void set_target_time(double target_time) override
Definition: transport.cc:670
int * get_row_4_el() const
Definition: mesh.h:159
double ** tm_diag
Definition: transport.h:286
virtual void output_data() override
Write computed fields.
Definition: transport.cc:883
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:149
Distributed sparse graphs, partitioning.
Vec vcfl_source_
Parallel vector for source term contribution to CFL condition.
Definition: transport.h:276
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:90
int * get_row_4_el() override
Return global array of order of elements within parallel vector.
Definition: transport.cc:856
FieldCommon & name(const string &name)
Definition: field_common.hh:86
#define END_TIMER(tag)
Ends a timer with specified tag.
arma::vec::fixed< spacedim > centre() const
Definition: accessors.hh:91
#define OLD_ASSERT_EQUAL(a, b)
Definition: global_defs.h:133
#define FOR_ELM_NEIGHS_VB(i, j)
Definition: elements.h:190
ElementAccessor< 3 > element_accessor() const
Gets ElementAccessor of this element.
Definition: elements.cc:148
double ** cumulative_corr
Definition: transport.h:303
bool changed() const
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:173
Class used for marking specified times at which some events occur.
Definition: time_marks.hh:36
Record type proxy class.
Definition: type_record.hh:177
void get_par_info(int *&el_4_loc, Distribution *&el_ds) override
Return array of indices of local elements and parallel distribution of elements.
Definition: transport.cc:846
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Class for representation SI units of Fields.
Definition: unit_si.hh:40
#define MPI_Barrier(comm)
Definition: mpi.h:531
double last_dt() const
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:240
Other possible transformation of coordinates:
Distribution * el_ds
Definition: transport.h:315
SideIter side(const unsigned int i) const
Definition: edges.h:31
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:47
void make_transport_partitioning()
Definition: transport.cc:152
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:274
TimeGovernor * time_
Definition: equation.hh:224
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:216
unsigned int lsize(int proc) const
get local size