Flow123d  release_2.2.0-48-gb04af7f
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, "Finite volume method, explicit in time, 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  "Specification of output fields and output times.")
72  .close();
73 }
74 
75 
77 {
78  ADD_FIELD(bc_conc, "Boundary condition for concentration of substances.", "0.0");
79  bc_conc.units( UnitSI().kg().m(-3) );
80  ADD_FIELD(init_conc, "Initial values for concentration of substances.", "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  .description("Concentration solution.");
88  output_fields += region_id.name("region_id")
91  .description("Region ids.");
92  output_fields += subdomain.name("subdomain")
95  .description("Subdomain ids of the domain decomposition.");
96 }
97 
98 
100 : ConcentrationTransportBase(init_mesh, in_rec),
101  is_mass_diag_changed(false),
102  input_rec(in_rec),
103  mh_dh(nullptr),
104  sources_corr(nullptr)
105 {
106  START_TIMER("ConvectionTransport");
107  this->eq_data_ = &data_;
108 
109  transport_matrix_time = -1.0; // or -infty
110  transport_bc_time = -1.0;
112  is_src_term_scaled = false;
113  is_bc_term_scaled = false;
114 }
115 
117 {
119 
121 
123  data_.set_input_list( input_rec.val<Input::Array>("input_fields") );
124  data_.set_mesh(*mesh_);
125 
129 
130  // register output vectors
137  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
138  {
139  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
140  auto output_field_ptr = out_conc[sbi].create_field<3, FieldValue<3>::Scalar>(n_substances());
141  data_.conc_mobile[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
142  }
143  //output_stream_->add_admissible_field_names(input_rec.val<Input::Array>("output_fields"));
144  //output_stream_->mark_output_times(*time_);
145 
147  //cout << "Transport." << endl;
148  //cout << time().marks();
149 
150  balance_->allocate(el_ds->lsize(), 1);
151 
152 }
153 
154 
155 //=============================================================================
156 // MAKE TRANSPORT
157 //=============================================================================
159 
160 // int * id_4_old = new int[mesh_->n_elements()];
161 // int i = 0;
162 // FOR_ELEMENTS(mesh_, ele) id_4_old[i++] = ele.index();
163 // mesh_->get_part()->id_maps(mesh_->n_elements(), id_4_old, el_ds, el_4_loc, row_4_el);
164 // delete[] id_4_old;
165  el_ds = mesh_->get_el_ds();
168 
169  // TODO: make output of partitioning is usefull but makes outputs different
170  // on different number of processors, which breaks tests.
171  //
172  // Possible solution:
173  // - have flag in ini file to turn this output ON
174  // - possibility to have different ref_output for different num of proc.
175  // - or do not test such kind of output
176  //
177  //FOR_ELEMENTS(mesh_, ele) {
178  // ele->pid=el_ds->get_proc(row_4_el[ele.index()]);
179  //}
180 
181 }
182 
183 
184 
186 {
187  unsigned int sbi;
188 
189  if (sources_corr) {
190  //Destroy mpi vectors at first
191  MatDestroy(&tm);
192  VecDestroy(&mass_diag);
193  VecDestroy(&vpmass_diag);
194  VecDestroy(&vcfl_flow_);
195  VecDestroy(&vcfl_source_);
196  delete cfl_flow_;
197  delete cfl_source_;
198 
199  for (sbi = 0; sbi < n_substances(); sbi++) {
200  // mpi vectors
201  VecDestroy(&(vconc[sbi]));
202  VecDestroy(&(vpconc[sbi]));
203  VecDestroy(&(bcvcorr[sbi]));
204  VecDestroy(&(vcumulative_corr[sbi]));
205  VecDestroy(&(v_tm_diag[sbi]));
206  VecDestroy(&(v_sources_corr[sbi]));
207 
208  // arrays of arrays
209  delete conc[sbi];
210  delete cumulative_corr[sbi];
211  delete tm_diag[sbi];
212  delete sources_corr[sbi];
213  }
214 
215  // arrays of mpi vectors
216  delete vconc;
217  delete vpconc;
218  delete bcvcorr;
219  delete vcumulative_corr;
220  delete v_tm_diag;
221  delete v_sources_corr;
222 
223  // arrays of arrays
224  delete conc;
225  delete cumulative_corr;
226  delete tm_diag;
227  delete sources_corr;
228  }
229 }
230 
231 
232 
233 
234 
236 {
237  FOR_ELEMENTS(mesh_, elem)
238  {
239  if (!el_ds->is_local(row_4_el[elem.index()])) continue;
240 
241  unsigned int index = row_4_el[elem.index()] - el_ds->begin();
242  ElementAccessor<3> ele_acc = mesh_->element_accessor(elem.index());
243 
244  for (unsigned int sbi=0; sbi<n_substances(); sbi++) // Optimize: SWAP LOOPS
245  conc[sbi][index] = data_.init_conc[sbi].value(elem->centre(), ele_acc);
246  }
247 
248 }
249 
250 //=============================================================================
251 // ALLOCATE OF TRANSPORT VARIABLES (ELEMENT & NODES)
252 //=============================================================================
254 
255  unsigned int i, sbi, n_subst;
256  n_subst = n_substances();
257 
258  sources_corr = new double*[n_subst];
259  tm_diag = new double*[n_subst];
260  cumulative_corr = new double*[n_subst];
261  for (sbi = 0; sbi < n_subst; sbi++) {
262  cumulative_corr[sbi] = new double[el_ds->lsize()];
263  sources_corr[sbi] = new double[el_ds->lsize()];
264  tm_diag[sbi] = new double[el_ds->lsize()];
265  }
266 
267  conc = new double*[n_subst];
268  out_conc.clear();
269  out_conc.resize(n_subst);
270  for (sbi = 0; sbi < n_subst; sbi++) {
271  conc[sbi] = new double[el_ds->lsize()];
272  out_conc[sbi].resize( el_ds->size() );
273  for (i = 0; i < el_ds->lsize(); i++) {
274  conc[sbi][i] = 0.0;
275  }
276  }
277 
278  cfl_flow_ = new double[el_ds->lsize()];
279  cfl_source_ = new double[el_ds->lsize()];
280 }
281 
282 //=============================================================================
283 // ALLOCATION OF TRANSPORT VECTORS (MPI)
284 //=============================================================================
286 
287  int sbi, n_subst, rank, np;
288  n_subst = n_substances();
289 
290  MPI_Barrier(PETSC_COMM_WORLD);
291  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
292  MPI_Comm_size(PETSC_COMM_WORLD, &np);
293 
294  vconc = new Vec[n_subst];
295  vpconc = new Vec[n_subst];
296  bcvcorr = new Vec[n_subst];
297  vcumulative_corr = new Vec[n_subst];
298  v_tm_diag = new Vec[n_subst];
299  v_sources_corr = new Vec[n_subst];
300 
301 
302  for (sbi = 0; sbi < n_subst; sbi++) {
303  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &bcvcorr[sbi]);
304  VecZeroEntries(bcvcorr[sbi]);
305  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(), conc[sbi],
306  &vconc[sbi]);
307 
308  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpconc[sbi]);
309  VecZeroEntries(vconc[sbi]);
310  VecZeroEntries(vpconc[sbi]);
311 
312  // SOURCES
313  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
314  cumulative_corr[sbi],&vcumulative_corr[sbi]);
315 
316  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
317  sources_corr[sbi],&v_sources_corr[sbi]);
318 
319  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
320  tm_diag[sbi],&v_tm_diag[sbi]);
321 
322  VecZeroEntries(vcumulative_corr[sbi]);
323  VecZeroEntries(out_conc[sbi].get_data_petsc());
324  }
325 
326 
327  MatCreateAIJ(PETSC_COMM_WORLD, el_ds->lsize(), el_ds->lsize(), mesh_->n_elements(),
328  mesh_->n_elements(), 16, PETSC_NULL, 4, PETSC_NULL, &tm);
329 
330  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &mass_diag);
331  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpmass_diag);
332 
333  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
335  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
337 }
338 
339 
341 {
342  START_TIMER ("set_boundary_conditions");
343 
345 
346  unsigned int sbi, loc_el, loc_b = 0;
347 
348  // Assembly bcvcorr vector
349  for(sbi=0; sbi < n_substances(); sbi++) VecZeroEntries(bcvcorr[sbi]);
350 
351  balance_->start_flux_assembly(subst_idx);
352 
353  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
354  elm = mesh_->element(el_4_loc[loc_el]);
355  if (elm->boundary_idx_ != NULL) {
356  unsigned int new_i = row_4_el[elm.index()];
357 
358  FOR_ELEMENT_SIDES(elm,si) {
359  Boundary *b = elm->side(si)->cond();
360  if (b != NULL) {
361  double flux = mh_dh->side_flux( *(elm->side(si)) );
362  if (flux < 0.0) {
363  double aij = -(flux / elm->measure() );
364 
365  for (sbi=0; sbi<n_substances(); sbi++)
366  {
367  double value = data_.bc_conc[sbi].value( b->element()->centre(), b->element_accessor() );
368 
369  VecSetValue(bcvcorr[sbi], new_i, value * aij, ADD_VALUES);
370 
371  // CAUTION: It seems that PETSc possibly optimize allocated space during assembly.
372  // So we have to add also values that may be non-zero in future due to changing velocity field.
373  balance_->add_flux_matrix_values(subst_idx[sbi], loc_b, {row_4_el[el_4_loc[loc_el]]}, {0.});
374  balance_->add_flux_vec_value(subst_idx[sbi], loc_b, flux*value);
375  }
376  } else {
377  for (sbi=0; sbi<n_substances(); sbi++)
378  VecSetValue(bcvcorr[sbi], new_i, 0, ADD_VALUES);
379 
380  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
381  {
382  balance_->add_flux_matrix_values(subst_idx[sbi], loc_b, {row_4_el[el_4_loc[loc_el]]}, {flux});
383  balance_->add_flux_vec_value(subst_idx[sbi], loc_b, 0);
384  }
385  }
386  ++loc_b;
387  }
388  }
389 
390  }
391  }
392 
393  balance_->finish_flux_assembly(subst_idx);
394 
395  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyBegin(bcvcorr[sbi]);
396  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyEnd(bcvcorr[sbi]);
397 
398  // we are calling set_boundary_conditions() after next_time() and
399  // we are using data from t() before, so we need to set corresponding bc time
401 }
402 
403 
404 //=============================================================================
405 // COMPUTE SOURCES
406 //=============================================================================
408 
409  //temporary variables
410  unsigned int loc_el, sbi;
411  double csection, source, diag;
412 
413  Element *ele;
414  ElementAccessor<3> ele_acc;
415  arma::vec3 p;
416 
417  //TODO: would it be possible to check the change in data for chosen substance? (may be in multifields?)
418 
419  //checking if the data were changed
420  if( (data_.sources_density.changed() )
421  || (data_.sources_conc.changed() )
422  || (data_.sources_sigma.changed() )
423  || (data_.cross_section.changed()))
424  {
425  START_TIMER("sources_reinit");
426  balance_->start_source_assembly(subst_idx);
427 
428  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++)
429  {
430  ele = mesh_->element(el_4_loc[loc_el]);
431  ele_acc = ele->element_accessor();
432  p = ele_acc.centre();
433  csection = data_.cross_section.value(p, ele_acc);
434 
435  // read for all substances
436  double max_cfl=0;
437  for (sbi = 0; sbi < n_substances(); sbi++)
438  {
439  double src_sigma = data_.sources_sigma[sbi].value(p, ele_acc);
440 
441  source = csection * (data_.sources_density[sbi].value(p, ele_acc) + src_sigma * data_.sources_conc[sbi].value(p, ele_acc));
442  // addition to RHS
443  sources_corr[sbi][loc_el] = source;
444  // addition to diagonal of the transport matrix
445  diag = src_sigma * csection;
446  tm_diag[sbi][loc_el] = - diag;
447 
448  // compute maximal cfl condition over all substances
449  max_cfl = std::max(max_cfl, fabs(diag));
450 
451  balance_->add_source_matrix_values(sbi, ele_acc.region().bulk_idx(), {row_4_el[el_4_loc[loc_el]]},
452  {- src_sigma * ele->measure() * csection});
453  balance_->add_source_vec_values(sbi, ele_acc.region().bulk_idx(), {row_4_el[el_4_loc[loc_el]]},
454  {source * ele->measure()});
455  }
456 
457  cfl_source_[loc_el] = max_cfl;
458  }
459 
460  balance_->finish_source_assembly(subst_idx);
461 
462  END_TIMER("sources_reinit");
463  }
464 }
465 
466 
467 
469 {
471 
474  std::stringstream ss; // print warning message with table of uninitialized fields
475  if ( FieldCommon::print_message_table(ss, "convection transport") ) {
476  WarningOut() << ss.str();
477  }
478 
481 
482  START_TIMER("Convection balance zero time step");
483 
487 
488  // write initial condition
489  output_data();
490 }
491 
492 
494 {
495  OLD_ASSERT(mh_dh, "Null MH object.\n" );
496  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
497 
498  START_TIMER("data reinit");
499 
500  bool cfl_changed = false;
501 
502  // if FLOW or DATA changed ---------------------> recompute transport matrix
504  {
507  cfl_changed = true;
508  DebugOut() << "CFL changed - flow.\n";
509  }
510 
512  {
514  cfl_changed = true;
515  DebugOut() << "CFL changed - mass matrix.\n";
516  }
517 
518  // if DATA changed ---------------------> recompute concentration sources (rhs and matrix diagonal)
521  {
523  is_src_term_scaled = false;
524  cfl_changed = true;
525  DebugOut() << "CFL changed - source.\n";
526  }
527 
528  // now resolve the CFL condition
529  if(cfl_changed)
530  {
531  // find maximum of sum of contribution from flow and sources: MAX(vcfl_flow_ + vcfl_source_)
532  Vec cfl;
533  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(),PETSC_DETERMINE, &cfl);
534  VecWAXPY(cfl, 1.0, vcfl_flow_, vcfl_source_);
535  VecMaxPointwiseDivide(cfl,mass_diag, &cfl_max_step);
536  // get a reciprocal value as a time constraint
538  DebugOut().fmt("CFL constraint (transport): {}\n", cfl_max_step);
539  }
540 
541  // although it does not influence CFL, compute BC so the full system is assembled
543  || data_.porosity.changed()
545  || data_.bc_conc.changed() )
546  {
548  is_bc_term_scaled = false;
549  }
550 
551  END_TIMER("data reinit");
552 
553  // return time constraint
554  time_constraint = cfl_max_step;
555  return cfl_changed;
556 }
557 
559 
560  START_TIMER("convection-one step");
561 
562  // proceed to next time - which we are about to compute
563  // explicit scheme looks one step back and uses data from previous time
564  // (data time set previously in assess_time_constraint())
565  time_->next_time();
566 
567  double dt_new = time_->dt(), // current time step we are about to compute
568  dt_scaled = dt_new / time_->last_dt(); // scaling ratio to previous time step
569 
570  START_TIMER("time step rescaling");
571 
572  // if FLOW or DATA or BC or DT changed ---------------------> rescale boundary condition
574  {
575  DebugOut() << "BC - rescale dt.\n";
576  //choose between fresh scaling with new dt or rescaling to a new dt
577  double dt = (!is_bc_term_scaled) ? dt_new : dt_scaled;
578  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
579  VecScale(bcvcorr[sbi], dt);
580  is_bc_term_scaled = true;
581  }
582 
583 
584  // if DATA or TIME STEP changed -----------------------> rescale source term
586  DebugOut() << "SRC - rescale dt.\n";
587  //choose between fresh scaling with new dt or rescaling to a new dt
588  double dt = (!is_src_term_scaled) ? dt_new : dt_scaled;
589  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
590  {
591  VecScale(v_sources_corr[sbi], dt);
592  VecScale(v_tm_diag[sbi], dt);
593  }
594  is_src_term_scaled = true;
595  }
596 
597  // if DATA or TIME STEP changed -----------------------> rescale transport matrix
599  DebugOut() << "TM - rescale dt.\n";
600  //choose between fresh scaling with new dt or rescaling to a new dt
601  double dt = (!is_convection_matrix_scaled) ? dt_new : dt_scaled;
602 
603  MatScale(tm, dt);
605  }
606 
607  END_TIMER("time step rescaling");
608 
609 
610  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
612  {
613  VecCopy(mass_diag, vpmass_diag);
615  } else is_mass_diag_changed = false;
616 
617 
618  // Compute new concentrations for every substance.
619 
620  for (unsigned int sbi = 0; sbi < n_substances(); sbi++) {
621  // one step in MOBILE phase
622  START_TIMER("mat mult");
623 
624  // tm_diag is a diagonal part of transport matrix, which depends on substance data (sources_sigma)
625  // Wwe need keep transport matrix independent of substance, therefore we keep this diagonal part
626  // separately in a vector tm_diag.
627  // Computation: first, we compute this diagonal addition D*pconc and save it temporaly into RHS
628 
629  // RHS = D*pconc, where D is diagonal matrix represented by a vector
630  VecPointwiseMult(vcumulative_corr[sbi], v_tm_diag[sbi], vconc[sbi]); //w = x.*y
631 
632  // Then we add boundary terms ans other source terms into RHS.
633  // RHS = 1.0 * bcvcorr + 1.0 * v_sources_corr + 1.0 * rhs
634  VecAXPBYPCZ(vcumulative_corr[sbi], 1.0, 1.0, 1.0, bcvcorr[sbi], v_sources_corr[sbi]); //z = ax + by + cz
635 
636  // Then we set the new previous concentration.
637  VecCopy(vconc[sbi], vpconc[sbi]); // pconc = conc
638  // And finally proceed with transport matrix multiplication.
639  if (is_mass_diag_changed) {
640  VecPointwiseMult(vconc[sbi], vconc[sbi], vpmass_diag); // vconc*=vpmass_diag
641  MatMultAdd(tm, vpconc[sbi], vconc[sbi], vconc[sbi]); // vconc+=tm*vpconc
642  VecAXPY(vconc[sbi], 1, vcumulative_corr[sbi]); // vconc+=vcumulative_corr
643  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
644  } else {
645  MatMultAdd(tm, vpconc[sbi], vcumulative_corr[sbi], vconc[sbi]); // vconc =tm*vpconc+vcumulative_corr
646  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
647  VecAXPY(vconc[sbi], 1, vpconc[sbi]); // vconc+=vpconc
648  }
649 
650  END_TIMER("mat mult");
651  }
652 
653  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
654  balance_->calculate_cumulative(sbi, vpconc[sbi]);
655 
656  END_TIMER("convection-one step");
657 }
658 
659 
660 void ConvectionTransport::set_target_time(double target_time)
661 {
662 
663  //sets target_mark_type (it is fixed) to be met in next_time()
664  time_->marks().add(TimeMark(target_time, target_mark_type));
665 
666  // This is done every time TOS calls update_solution.
667  // If CFL condition is changed, time fixation will change later from TOS.
668 
669  // Set the same constraint as was set last time.
670 
671  // TODO: fix this hack, remove this method completely, leaving just the first line at the calling point
672  // in TransportOperatorSplitting::update_solution()
673  // doing this directly leads to choose of large time step violationg CFL condition
674  if (cfl_max_step > time_->dt()*1e-10)
675  time_->set_upper_constraint(cfl_max_step, "CFL condition used from previous step.");
676 
677  // fixing convection time governor till next target_mark_type (got from TOS or other)
678  // may have marks for data changes
680 }
681 
682 
684 {
686 
687  VecZeroEntries(mass_diag);
688 
689  balance_->start_mass_assembly(subst_idx);
690 
691  for (unsigned int loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
692  elm = mesh_->element(el_4_loc[loc_el]);
693 
694  double csection = data_.cross_section.value(elm->centre(), elm->element_accessor());
695  //double por_m = data_.porosity.value(elm->centre(), elm->element_accessor());
696  double por_m = data_.water_content.value(elm->centre(), elm->element_accessor());
697 
698  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
699  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()} );
700 
701  VecSetValue(mass_diag, row_4_el[elm.index()], csection*por_m, INSERT_VALUES);
702  }
703 
704  balance_->finish_mass_assembly(subst_idx);
705 
706  VecAssemblyBegin(mass_diag);
707  VecAssemblyEnd(mass_diag);
708 
709  is_mass_diag_changed = true;
710 }
711 
712 
713 //=============================================================================
714 // CREATE TRANSPORT MATRIX
715 //=============================================================================
717 
718  START_TIMER("convection_matrix_assembly");
719 
722  struct Edge *edg;
723  unsigned int n;
724  int s, j, np, rank, new_j, new_i;
725  double aij, aii;
726 
727  MatZeroEntries(tm);
728 
729  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
730  MPI_Comm_size(PETSC_COMM_WORLD, &np);
731 
732  double flux, flux2, edg_flux;
733 
734  aii = 0.0;
735 
736  for (unsigned int loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
737  elm = mesh_->element(el_4_loc[loc_el]);
738  new_i = row_4_el[elm.index()];
739 
740  FOR_ELEMENT_SIDES(elm,si) {
741  // same dim
742  flux = mh_dh->side_flux( *(elm->side(si)) );
743  if (elm->side(si)->cond() == NULL) {
744  edg = elm->side(si)->edge();
745  edg_flux = 0;
746  for( int s=0; s < edg->n_sides; s++) {
747  flux2 = mh_dh->side_flux( *(edg->side(s)) );
748  if ( flux2 > 0) edg_flux+= flux2;
749  }
750  FOR_EDGE_SIDES(edg,s)
751  // this test should also eliminate sides facing to lower dim. elements in comp. neighboring
752  // These edges on these sides should have just one side
753  if (edg->side(s) != elm->side(si)) {
754  j = ELEMENT_FULL_ITER(mesh_, edg->side(s)->element()).index();
755  new_j = row_4_el[j];
756 
757  flux2 = mh_dh->side_flux( *(edg->side(s)));
758  if ( flux2 > 0.0 && flux <0.0)
759  aij = -(flux * flux2 / ( edg_flux * elm->measure() ) );
760  else aij =0;
761  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
762  }
763  }
764  if (flux > 0.0)
765  aii -= (flux / elm->measure() );
766  } // end same dim //ELEMENT_SIDES
767 
768  FOR_ELM_NEIGHS_VB(elm,n) // comp model
769  {
770  el2 = ELEMENT_FULL_ITER(mesh_, elm->neigh_vb[n]->side()->element() ); // higher dim. el.
771  OLD_ASSERT( el2 != elm, "Elm. same\n");
772  new_j = row_4_el[el2.index()];
773  flux = mh_dh->side_flux( *(elm->neigh_vb[n]->side()) );
774 
775  // volume source - out-flow from higher dimension
776  if (flux > 0.0) aij = flux / elm->measure();
777  else aij=0;
778  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
779  // out flow from higher dim. already accounted
780 
781  // volume drain - in-flow to higher dimension
782  if (flux < 0.0) {
783  aii -= (-flux) / elm->measure(); // diagonal drain
784  aij = (-flux) / el2->measure();
785  } else aij=0;
786  MatSetValue(tm, new_j, new_i, aij, INSERT_VALUES);
787  }
788 
789  MatSetValue(tm, new_i, new_i, aii, INSERT_VALUES);
790 
791  cfl_flow_[loc_el] = fabs(aii);
792  aii = 0.0;
793  } // END ELEMENTS
794 
795  MatAssemblyBegin(tm, MAT_FINAL_ASSEMBLY);
796  MatAssemblyEnd(tm, MAT_FINAL_ASSEMBLY);
797 
799  END_TIMER("convection_matrix_assembly");
800 
802 }
803 
804 
805 
806 
807 
808 //=============================================================================
809 // OUTPUT VECTOR GATHER
810 //=============================================================================
812 
813  unsigned int sbi;
814  IS is;
815 
816  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el, PETSC_COPY_VALUES, &is); //WithArray
817  VecScatterCreate(vconc[0], is, out_conc[0].get_data_petsc(), PETSC_NULL, &vconc_out_scatter);
818  for (sbi = 0; sbi < n_substances(); sbi++) {
819  VecScatterBegin(vconc_out_scatter, vconc[sbi], out_conc[sbi].get_data_petsc(), INSERT_VALUES, SCATTER_FORWARD);
820  VecScatterEnd(vconc_out_scatter, vconc[sbi], out_conc[sbi].get_data_petsc(), INSERT_VALUES, SCATTER_FORWARD);
821  }
822  VecScatterDestroy(&(vconc_out_scatter));
823  ISDestroy(&(is));
824 }
825 
826 
828  return conc;
829 }
830 
831 void ConvectionTransport::get_par_info(int * &el_4_loc_out, Distribution * &el_distribution_out){
832  el_4_loc_out = this->el_4_loc;
833  el_distribution_out = this->el_ds;
834  return;
835 }
836 
837 //int *ConvectionTransport::get_el_4_loc(){
838 // return el_4_loc;
839 //}
840 
842  return row_4_el;
843 }
844 
845 
846 
848 
850  //if ( data_.output_fields.is_field_output_time(data_.conc_mobile, time().step()) ) {
852  //}
853 
854  data_.output_fields.output(time().step());
855 
856  START_TIMER("TOS-balance");
857  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
858  balance_->calculate_instant(sbi, vconc[sbi]);
859  balance_->output();
860  END_TIMER("TOS-balance");
861 }
862 
863 void ConvectionTransport::set_balance_object(std::shared_ptr<Balance> balance)
864 {
865  balance_ = balance;
866  subst_idx = balance_->add_quantities(substances_.names());
867 }
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
FieldSet * eq_data_
Definition: equation.hh:232
static auto subdomain(Mesh &mesh) -> IndexField
std::vector< VectorSeqDouble > out_conc
Definition: transport.h:294
bool is_mass_diag_changed
Flag indicates that porosity or cross_section changed during last time.
Definition: transport.h:255
double measure() const
Definition: elements.cc:89
double time_changed() const
void update_solution() override
Definition: transport.cc:558
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
double end_time() const
End time.
double transport_matrix_time
Definition: transport.h:280
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
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:61
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:426
void alloc_transport_vectors()
Definition: transport.cc:253
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:189
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:90
void create_mass_matrix()
Definition: transport.cc:683
void next_time()
Proceed to the next time according to current estimated time step.
void initialize() override
Definition: transport.cc:116
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:299
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:262
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
void set_initial_condition()
Definition: transport.cc:235
double transport_bc_time
Time of the last update of the boundary condition terms.
Definition: transport.h:281
double * cfl_flow_
Definition: transport.h:267
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:340
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
#define ELEMENT_FULL_ITER(_mesh_, i)
Definition: mesh.h:434
Definition: mesh.h:97
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:94
void initialize(std::shared_ptr< OutputTime > stream, Input::Record in_rec, const TimeGovernor &tg)
double ** sources_corr
Definition: transport.h:258
#define FOR_ELEMENT_SIDES(i, j)
Definition: elements.h:188
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:180
int n_sides
Definition: edges.h:36
Definition: edges.h:26
const RegionDB & region_db() const
Definition: mesh.h:155
void set_balance_object(std::shared_ptr< Balance > balance) override
Definition: transport.cc:863
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:669
double ** conc
Concentrations for phase, substance, element.
Definition: transport.h:286
const TimeStep & step(int index=-1) const
virtual ~ConvectionTransport()
Definition: transport.cc:185
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:468
static TimeMarks & marks()
void output_vector_gather()
Definition: transport.cc:811
Basic time management class.
Class for declaration of inputs sequences.
Definition: type_base.hh:345
const Input::Record input_rec
Record with input specification.
Definition: transport.h:297
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:314
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:195
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void setup_components()
unsigned int n_elements() const
Definition: mesh.h:141
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:292
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:99
static const int registrar
Registrar of class to factory.
Definition: transport.h:239
bool is_changed_dt() const
#define START_TIMER(tag)
Starts a timer with specified tag.
SubstanceList substances_
Transported substances.
Definition: transport.h:307
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:362
double ** get_concentration_matrix() override
Getter for array of concentrations per element.
Definition: transport.cc:827
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:490
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
Definition: transport.h:317
Distribution * get_el_ds() const
Definition: mesh.h:174
#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:284
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:81
Field< 3, FieldValue< 3 >::Integer > region_id
Definition: transport.h:88
void create_transport_matrix_mpi()
Definition: transport.cc:716
Support classes for parallel programing.
void alloc_transport_structs_mpi()
Definition: transport.cc:285
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
FieldCommon & description(const string &description)
bool is_convection_matrix_scaled
Definition: transport.h:252
FLOW123D_FORCE_LINK_IN_CHILD(convectionTransport)
double * cfl_source_
Definition: transport.h:267
#define ELEMENT_FULL_ITER_NULL(_mesh_)
Definition: mesh.h:440
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:407
bool evaluate_time_constraint(double &time_constraint) override
Definition: transport.cc:493
VecScatter vconc_out_scatter
Definition: transport.h:270
double dt() const
Vec vcfl_flow_
Parallel vector for flow contribution to CFL condition.
Definition: transport.h:265
arma::vec3 centre() const
Definition: elements.cc:120
void set_target_time(double target_time) override
Definition: transport.cc:660
int * get_row_4_el() const
Definition: mesh.h:177
double ** tm_diag
Definition: transport.h:275
virtual void output_data() override
Write computed fields.
Definition: transport.cc:847
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:265
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:841
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:236
FieldCommon & name(const string &name)
Definition: field_common.hh:97
#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:189
ElementAccessor< 3 > element_accessor() const
Gets ElementAccessor of this element.
Definition: elements.cc:148
double ** cumulative_corr
Definition: transport.h:292
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:182
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:831
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
Field< 3, FieldValue< 3 >::Integer > subdomain
Definition: transport.h:89
#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:242
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:91
Other possible transformation of coordinates:
Distribution * el_ds
Definition: transport.h:304
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:158
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:263
TimeGovernor * time_
Definition: equation.hh:224
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:228
unsigned int lsize(int proc) const
get local size