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