Flow123d  release_3.0.0-1152-gdb4be9b
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/side_impl.hh"
25 #include "mesh/long_idx.hh"
26 #include "mesh/mesh.h"
27 #include "mesh/partitioning.hh"
28 #include "mesh/accessors.hh"
29 #include "mesh/range_wrapper.hh"
30 #include "mesh/neighbours.h"
31 #include "transport/transport.h"
32 
33 #include "la/distribution.hh"
34 
35 #include "la/sparse_graph.hh"
36 #include <iostream>
37 #include <iomanip>
38 #include <string>
39 
40 #include "io/output_time.hh"
41 #include "tools/time_governor.hh"
42 #include "coupling/balance.hh"
43 #include "input/accessors.hh"
44 #include "input/input_type.hh"
45 
47 #include "fields/field_values.hh"
48 #include "fields/field_fe.hh"
50 #include "fields/generic_field.hh"
51 
52 #include "reaction/isotherm.hh" // SorptionType enum
53 #include "flow/mh_dofhandler.hh"
54 
55 
56 FLOW123D_FORCE_LINK_IN_CHILD(convectionTransport);
57 
58 
59 namespace IT = Input::Type;
60 
61 const string _equation_name = "Solute_Advection_FV";
62 
64  Input::register_class< ConvectionTransport, Mesh &, const Input::Record >(_equation_name) +
66 
68 {
69  return IT::Record(_equation_name, "Finite volume method, explicit in time, for advection only solute transport.")
71  .declare_key("input_fields", IT::Array(
72  EqData().make_field_descriptor_type(_equation_name)),
74  "")
75  .declare_key("output",
76  EqData().output_fields.make_output_type(_equation_name, ""),
77  IT::Default("{ \"fields\": [ \"conc\" ] }"),
78  "Specification of output fields and output times.")
79  .close();
80 }
81 
82 
84 {
85  *this += bc_conc.name("bc_conc")
86  .description("Boundary condition for concentration of substances.")
87  .input_default("0.0")
88  .units( UnitSI().kg().m(-3) );
89 
90  *this += init_conc.name("init_conc")
91  .description("Initial values for concentration of substances.")
92  .input_default("0.0")
93  .units( UnitSI().kg().m(-3) );
94 
95  output_fields += *this;
96  output_fields += conc_mobile.name("conc")
97  .units( UnitSI().kg().m(-3) )
99  .description("Concentration solution.");
100  output_fields += region_id.name("region_id")
103  .description("Region ids.");
104  output_fields += subdomain.name("subdomain")
107  .description("Subdomain ids of the domain decomposition.");
108 }
109 
110 
112 : ConcentrationTransportBase(init_mesh, in_rec),
113  is_mass_diag_changed(false),
114  input_rec(in_rec),
115  mh_dh(nullptr),
116  sources_corr(nullptr)
117 {
118  START_TIMER("ConvectionTransport");
119  this->eq_data_ = &data_;
120 
121  transport_matrix_time = -1.0; // or -infty
122  transport_bc_time = -1.0;
124  is_src_term_scaled = false;
125  is_bc_term_scaled = false;
126 
127  //initialization of DOF handler
128  static FE_P_disc<0> fe0(0);
129  static FE_P_disc<1> fe1(0);
130  static FE_P_disc<2> fe2(0);
131  static FE_P_disc<3> fe3(0);
132  dh_ = make_shared<DOFHandlerMultiDim>(init_mesh);
133  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>( &init_mesh, &fe0, &fe1, &fe2, &fe3);
134  dh_->distribute_dofs(ds);
135 
136 }
137 
139 {
141 
143 
145  data_.set_input_list( input_rec.val<Input::Array>("input_fields"), *time_ );
146  data_.set_mesh(*mesh_);
147 
151 
152  // register output vectors
159  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
160  {
161  // create shared pointer to a FieldFE and push this Field to output_field on all regions
162  output_field_ptr[sbi] = create_field<3, FieldValue<3>::Scalar>(out_conc[sbi], *mesh_, 1);
163  data_.conc_mobile[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr[sbi], 0);
164  }
165  //output_stream_->add_admissible_field_names(input_rec.val<Input::Array>("output_fields"));
166  //output_stream_->mark_output_times(*time_);
167 
169  //cout << "Transport." << endl;
170  //cout << time().marks();
171 
172  balance_->allocate(el_ds->lsize(), 1);
173 
174 }
175 
176 
177 //=============================================================================
178 // MAKE TRANSPORT
179 //=============================================================================
181 
182 // int * id_4_old = new int[mesh_->n_elements()];
183 // int i = 0;
184 // for (auto ele : mesh_->elements_range()) id_4_old[i++] = ele.index();
185 // mesh_->get_part()->id_maps(mesh_->n_elements(), id_4_old, el_ds, el_4_loc, row_4_el);
186 // delete[] id_4_old;
187  el_ds = mesh_->get_el_ds();
190 
191  // TODO: make output of partitioning is usefull but makes outputs different
192  // on different number of processors, which breaks tests.
193  //
194  // Possible solution:
195  // - have flag in ini file to turn this output ON
196  // - possibility to have different ref_output for different num of proc.
197  // - or do not test such kind of output
198  //
199  //for (auto ele : mesh_->elements_range()) {
200  // ele->pid()=el_ds->get_proc(row_4_el[ele.index()]);
201  //}
202 
203 }
204 
205 
206 
208 {
209  unsigned int sbi;
210 
211  if (sources_corr) {
212  //Destroy mpi vectors at first
213  chkerr(MatDestroy(&tm));
214  chkerr(VecDestroy(&mass_diag));
215  chkerr(VecDestroy(&vpmass_diag));
216  chkerr(VecDestroy(&vcfl_flow_));
217  chkerr(VecDestroy(&vcfl_source_));
218  delete cfl_flow_;
219  delete cfl_source_;
220 
221  for (sbi = 0; sbi < n_substances(); sbi++) {
222  // mpi vectors
223  chkerr(VecDestroy(&vconc[sbi]));
224  chkerr(VecDestroy(&vpconc[sbi]));
225  chkerr(VecDestroy(&bcvcorr[sbi]));
226  chkerr(VecDestroy(&vcumulative_corr[sbi]));
227  chkerr(VecDestroy(&v_tm_diag[sbi]));
228  chkerr(VecDestroy(&v_sources_corr[sbi]));
229 
230  // arrays of arrays
231  delete conc[sbi];
232  delete cumulative_corr[sbi];
233  delete tm_diag[sbi];
234  delete sources_corr[sbi];
235  }
236 
237  // arrays of mpi vectors
238  delete vconc;
239  delete vpconc;
240  delete bcvcorr;
241  delete vcumulative_corr;
242  delete v_tm_diag;
243  delete v_sources_corr;
244 
245  // arrays of arrays
246  delete conc;
247  delete cumulative_corr;
248  delete tm_diag;
249  delete sources_corr;
250  }
251 }
252 
253 
254 
255 
256 
258 {
259  for (auto elem : mesh_->elements_range()) {
260  if (!el_ds->is_local(row_4_el[ elem.idx() ])) continue;
261 
262  LongIdx index = row_4_el[ elem.idx() ] - el_ds->begin();
263  ElementAccessor<3> ele_acc = mesh_->element_accessor( elem.idx() );
264 
265  for (unsigned int sbi=0; sbi<n_substances(); sbi++) // Optimize: SWAP LOOPS
266  conc[sbi][index] = data_.init_conc[sbi].value(elem.centre(), ele_acc);
267  }
268 
269 }
270 
271 //=============================================================================
272 // ALLOCATE OF TRANSPORT VARIABLES (ELEMENT & NODES)
273 //=============================================================================
275 
276  unsigned int i, sbi, n_subst;
277  n_subst = n_substances();
278 
279  sources_corr = new double*[n_subst];
280  tm_diag = new double*[n_subst];
281  cumulative_corr = new double*[n_subst];
282  for (sbi = 0; sbi < n_subst; sbi++) {
283  cumulative_corr[sbi] = new double[el_ds->lsize()];
284  sources_corr[sbi] = new double[el_ds->lsize()];
285  tm_diag[sbi] = new double[el_ds->lsize()];
286  }
287 
288  conc = new double*[n_subst];
289  out_conc.clear();
290  out_conc.resize(n_subst);
291  output_field_ptr.clear();
292  output_field_ptr.resize(n_subst);
293  for (sbi = 0; sbi < n_subst; sbi++) {
294  conc[sbi] = new double[el_ds->lsize()];
295  out_conc[sbi].resize( el_ds->size() );
296  for (i = 0; i < el_ds->lsize(); i++) {
297  conc[sbi][i] = 0.0;
298  }
299  }
300 
301  cfl_flow_ = new double[el_ds->lsize()];
302  cfl_source_ = new double[el_ds->lsize()];
303 }
304 
305 //=============================================================================
306 // ALLOCATION OF TRANSPORT VECTORS (MPI)
307 //=============================================================================
309 
310  int sbi, n_subst;
311  n_subst = n_substances();
312 
313  MPI_Barrier(PETSC_COMM_WORLD);
314 
315  vconc = new Vec[n_subst];
316  vpconc = new Vec[n_subst];
317  bcvcorr = new Vec[n_subst];
318  vcumulative_corr = new Vec[n_subst];
319  v_tm_diag = new Vec[n_subst];
320  v_sources_corr = new Vec[n_subst];
321 
322 
323  for (sbi = 0; sbi < n_subst; sbi++) {
324  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &bcvcorr[sbi]);
325  VecZeroEntries(bcvcorr[sbi]);
326  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(), conc[sbi],
327  &vconc[sbi]);
328 
329  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpconc[sbi]);
330  VecZeroEntries(vconc[sbi]);
331  VecZeroEntries(vpconc[sbi]);
332 
333  // SOURCES
334  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
335  cumulative_corr[sbi],&vcumulative_corr[sbi]);
336 
337  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
338  sources_corr[sbi],&v_sources_corr[sbi]);
339 
340  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
341  tm_diag[sbi],&v_tm_diag[sbi]);
342 
343  VecZeroEntries(vcumulative_corr[sbi]);
344  VecZeroEntries(out_conc[sbi].petsc_vec());
345  }
346 
347 
348  MatCreateAIJ(PETSC_COMM_WORLD, el_ds->lsize(), el_ds->lsize(), mesh_->n_elements(),
349  mesh_->n_elements(), 16, PETSC_NULL, 4, PETSC_NULL, &tm);
350 
351  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &mass_diag);
352  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpmass_diag);
353 
354  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
356  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
358 }
359 
360 
362 {
363  START_TIMER ("set_boundary_conditions");
364 
365  ElementAccessor<3> elm;
366 
367  unsigned int sbi, loc_el;
368 
369  // Assembly bcvcorr vector
370  for(sbi=0; sbi < n_substances(); sbi++) VecZeroEntries(bcvcorr[sbi]);
371 
372  balance_->start_flux_assembly(subst_idx);
373 
374  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
375  elm = mesh_->element_accessor( el_4_loc[loc_el] );
376  if (elm->boundary_idx_ != NULL) {
377  LongIdx new_i = row_4_el[ elm.idx() ];
378 
379  for (unsigned int si=0; si<elm->n_sides(); si++) {
380  Boundary *b = elm.side(si)->cond();
381  if (b != NULL) {
382  double flux = mh_dh->side_flux( *(elm.side(si)) );
383  if (flux < 0.0) {
384  double aij = -(flux / elm.measure() );
385 
386  for (sbi=0; sbi<n_substances(); sbi++)
387  {
388  double value = data_.bc_conc[sbi].value( b->element_accessor().centre(), b->element_accessor() );
389 
390  VecSetValue(bcvcorr[sbi], new_i, value * aij, ADD_VALUES);
391 
392  // CAUTION: It seems that PETSc possibly optimize allocated space during assembly.
393  // So we have to add also values that may be non-zero in future due to changing velocity field.
394  balance_->add_flux_matrix_values(subst_idx[sbi], elm.side(si), {row_4_el[el_4_loc[loc_el]]}, {0.});
395  balance_->add_flux_vec_value(subst_idx[sbi], elm.side(si), flux*value);
396  }
397  } else {
398  for (sbi=0; sbi<n_substances(); sbi++)
399  VecSetValue(bcvcorr[sbi], new_i, 0, ADD_VALUES);
400 
401  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
402  {
403  balance_->add_flux_matrix_values(subst_idx[sbi], elm.side(si), {row_4_el[el_4_loc[loc_el]]}, {flux});
404  balance_->add_flux_vec_value(subst_idx[sbi], elm.side(si), 0);
405  }
406  }
407  }
408  }
409 
410  }
411  }
412 
413  balance_->finish_flux_assembly(subst_idx);
414 
415  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyBegin(bcvcorr[sbi]);
416  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyEnd(bcvcorr[sbi]);
417 
418  // we are calling set_boundary_conditions() after next_time() and
419  // we are using data from t() before, so we need to set corresponding bc time
421 }
422 
423 
424 //=============================================================================
425 // COMPUTE SOURCES
426 //=============================================================================
428 
429  //temporary variables
430  unsigned int loc_el, sbi;
431  double csection, source, diag;
432 
433  ElementAccessor<3> ele_acc;
434  arma::vec3 p;
435 
436  //TODO: would it be possible to check the change in data for chosen substance? (may be in multifields?)
437 
438  //checking if the data were changed
439  if( (data_.sources_density.changed() )
440  || (data_.sources_conc.changed() )
441  || (data_.sources_sigma.changed() )
442  || (data_.cross_section.changed()))
443  {
444  START_TIMER("sources_reinit");
445  balance_->start_source_assembly(subst_idx);
446 
447  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++)
448  {
449  ele_acc = mesh_->element_accessor( el_4_loc[loc_el] );
450  p = ele_acc.centre();
451  csection = data_.cross_section.value(p, ele_acc);
452 
453  // read for all substances
454  double max_cfl=0;
455  for (sbi = 0; sbi < n_substances(); sbi++)
456  {
457  double src_sigma = data_.sources_sigma[sbi].value(p, ele_acc);
458 
459  source = csection * (data_.sources_density[sbi].value(p, ele_acc) + src_sigma * data_.sources_conc[sbi].value(p, ele_acc));
460  // addition to RHS
461  sources_corr[sbi][loc_el] = source;
462  // addition to diagonal of the transport matrix
463  diag = src_sigma * csection;
464  tm_diag[sbi][loc_el] = - diag;
465 
466  // compute maximal cfl condition over all substances
467  max_cfl = std::max(max_cfl, fabs(diag));
468 
469  balance_->add_source_values(sbi, ele_acc.region().bulk_idx(), {int(loc_el)},
470  {- src_sigma * ele_acc.measure() * csection},
471  {source * ele_acc.measure()});
472  }
473 
474  cfl_source_[loc_el] = max_cfl;
475  }
476 
477  balance_->finish_source_assembly(subst_idx);
478 
479  END_TIMER("sources_reinit");
480  }
481 }
482 
483 
484 
486 {
488 
491  std::stringstream ss; // print warning message with table of uninitialized fields
492  if ( FieldCommon::print_message_table(ss, "convection transport") ) {
493  WarningOut() << ss.str();
494  }
495 
498 
499  START_TIMER("Convection balance zero time step");
500 
504 
505  // write initial condition
506  output_data();
507 }
508 
509 
511 {
512  OLD_ASSERT(mh_dh, "Null MH object.\n" );
513  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
514 
515  START_TIMER("data reinit");
516 
517  bool cfl_changed = false;
518 
519  // if FLOW or DATA changed ---------------------> recompute transport matrix
521  {
524  cfl_changed = true;
525  DebugOut() << "CFL changed - flow.\n";
526  }
527 
529  {
531  cfl_changed = true;
532  DebugOut() << "CFL changed - mass matrix.\n";
533  }
534 
535  // if DATA changed ---------------------> recompute concentration sources (rhs and matrix diagonal)
538  {
540  is_src_term_scaled = false;
541  cfl_changed = true;
542  DebugOut() << "CFL changed - source.\n";
543  }
544 
545  // now resolve the CFL condition
546  if(cfl_changed)
547  {
548  // find maximum of sum of contribution from flow and sources: MAX(vcfl_flow_ + vcfl_source_)
549  Vec cfl;
550  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(),PETSC_DETERMINE, &cfl);
551  VecWAXPY(cfl, 1.0, vcfl_flow_, vcfl_source_);
552  VecMaxPointwiseDivide(cfl,mass_diag, &cfl_max_step);
553  // get a reciprocal value as a time constraint
555  DebugOut().fmt("CFL constraint (transport): {}\n", cfl_max_step);
556  }
557 
558  // although it does not influence CFL, compute BC so the full system is assembled
560  || data_.porosity.changed()
562  || data_.bc_conc.changed() )
563  {
565  is_bc_term_scaled = false;
566  }
567 
568  END_TIMER("data reinit");
569 
570  // return time constraint
571  time_constraint = cfl_max_step;
572  return cfl_changed;
573 }
574 
576 
577  START_TIMER("convection-one step");
578 
579  // proceed to next time - which we are about to compute
580  // explicit scheme looks one step back and uses data from previous time
581  // (data time set previously in assess_time_constraint())
582  time_->next_time();
583 
584  double dt_new = time_->dt(), // current time step we are about to compute
585  dt_scaled = dt_new / time_->last_dt(); // scaling ratio to previous time step
586 
587  START_TIMER("time step rescaling");
588 
589  // if FLOW or DATA or BC or DT changed ---------------------> rescale boundary condition
591  {
592  DebugOut() << "BC - rescale dt.\n";
593  //choose between fresh scaling with new dt or rescaling to a new dt
594  double dt = (!is_bc_term_scaled) ? dt_new : dt_scaled;
595  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
596  VecScale(bcvcorr[sbi], dt);
597  is_bc_term_scaled = true;
598  }
599 
600 
601  // if DATA or TIME STEP changed -----------------------> rescale source term
603  DebugOut() << "SRC - rescale dt.\n";
604  //choose between fresh scaling with new dt or rescaling to a new dt
605  double dt = (!is_src_term_scaled) ? dt_new : dt_scaled;
606  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
607  {
608  VecScale(v_sources_corr[sbi], dt);
609  VecScale(v_tm_diag[sbi], dt);
610  }
611  is_src_term_scaled = true;
612  }
613 
614  // if DATA or TIME STEP changed -----------------------> rescale transport matrix
616  DebugOut() << "TM - rescale dt.\n";
617  //choose between fresh scaling with new dt or rescaling to a new dt
618  double dt = (!is_convection_matrix_scaled) ? dt_new : dt_scaled;
619 
620  MatScale(tm, dt);
622  }
623 
624  END_TIMER("time step rescaling");
625 
626 
627  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
629  {
630  VecCopy(mass_diag, vpmass_diag);
632  } else is_mass_diag_changed = false;
633 
634 
635  // Compute new concentrations for every substance.
636 
637  for (unsigned int sbi = 0; sbi < n_substances(); sbi++) {
638  // one step in MOBILE phase
639  START_TIMER("mat mult");
640 
641  // tm_diag is a diagonal part of transport matrix, which depends on substance data (sources_sigma)
642  // Wwe need keep transport matrix independent of substance, therefore we keep this diagonal part
643  // separately in a vector tm_diag.
644  // Computation: first, we compute this diagonal addition D*pconc and save it temporaly into RHS
645 
646  // RHS = D*pconc, where D is diagonal matrix represented by a vector
647  VecPointwiseMult(vcumulative_corr[sbi], v_tm_diag[sbi], vconc[sbi]); //w = x.*y
648 
649  // Then we add boundary terms ans other source terms into RHS.
650  // RHS = 1.0 * bcvcorr + 1.0 * v_sources_corr + 1.0 * rhs
651  VecAXPBYPCZ(vcumulative_corr[sbi], 1.0, 1.0, 1.0, bcvcorr[sbi], v_sources_corr[sbi]); //z = ax + by + cz
652 
653  // Then we set the new previous concentration.
654  VecCopy(vconc[sbi], vpconc[sbi]); // pconc = conc
655  // And finally proceed with transport matrix multiplication.
656  if (is_mass_diag_changed) {
657  VecPointwiseMult(vconc[sbi], vconc[sbi], vpmass_diag); // vconc*=vpmass_diag
658  MatMultAdd(tm, vpconc[sbi], vconc[sbi], vconc[sbi]); // vconc+=tm*vpconc
659  VecAXPY(vconc[sbi], 1, vcumulative_corr[sbi]); // vconc+=vcumulative_corr
660  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
661  } else {
662  MatMultAdd(tm, vpconc[sbi], vcumulative_corr[sbi], vconc[sbi]); // vconc =tm*vpconc+vcumulative_corr
663  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
664  VecAXPY(vconc[sbi], 1, vpconc[sbi]); // vconc+=vpconc
665  }
666 
667  END_TIMER("mat mult");
668  }
669 
670  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
671  balance_->calculate_cumulative(sbi, vpconc[sbi]);
672 
673  END_TIMER("convection-one step");
674 }
675 
676 
677 void ConvectionTransport::set_target_time(double target_time)
678 {
679 
680  //sets target_mark_type (it is fixed) to be met in next_time()
681  time_->marks().add(TimeMark(target_time, target_mark_type));
682 
683  // This is done every time TOS calls update_solution.
684  // If CFL condition is changed, time fixation will change later from TOS.
685 
686  // Set the same constraint as was set last time.
687 
688  // TODO: fix this hack, remove this method completely, leaving just the first line at the calling point
689  // in TransportOperatorSplitting::update_solution()
690  // doing this directly leads to choose of large time step violationg CFL condition
691  if (cfl_max_step > time_->dt()*1e-10)
692  time_->set_upper_constraint(cfl_max_step, "CFL condition used from previous step.");
693 
694  // fixing convection time governor till next target_mark_type (got from TOS or other)
695  // may have marks for data changes
697 }
698 
699 
701 {
702  ElementAccessor<3> elm;
703 
704  VecZeroEntries(mass_diag);
705 
706  balance_->start_mass_assembly(subst_idx);
707 
708  for (unsigned int loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
709  elm = mesh_->element_accessor( el_4_loc[loc_el] );
710 
711  double csection = data_.cross_section.value(elm.centre(), elm);
712  //double por_m = data_.porosity.value(elm.centre(), elm->element_accessor());
713  double por_m = data_.water_content.value(elm.centre(), elm);
714 
715  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
716  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()} );
717 
718  VecSetValue(mass_diag, row_4_el[ elm.idx() ], csection*por_m, INSERT_VALUES);
719  }
720 
721  balance_->finish_mass_assembly(subst_idx);
722 
723  VecAssemblyBegin(mass_diag);
724  VecAssemblyEnd(mass_diag);
725 
726  is_mass_diag_changed = true;
727 }
728 
729 
730 //=============================================================================
731 // CREATE TRANSPORT MATRIX
732 //=============================================================================
734 
735  START_TIMER("convection_matrix_assembly");
736 
737  ElementAccessor<3> el2;
738  ElementAccessor<3> elm;
739  const Edge *edg;
740  int j;
741  LongIdx new_j, new_i;
742  double aij, aii;
743 
744  MatZeroEntries(tm);
745 
746  double flux, flux2, edg_flux;
747 
748  aii = 0.0;
749 
750  unsigned int loc_el = 0;
751  for ( DHCellAccessor dh_cell : dh_->own_range() ) {
752  new_i = row_4_el[ dh_cell.elm_idx() ];
753  for( DHCellSide cell_side : dh_cell.side_range() ) {
754  flux = mh_dh->side_flux( cell_side.side());
755  if (cell_side.cond() == NULL) {
756  edg_flux = 0;
757  for( DHCellSide edge_side : cell_side.edge_sides() ) {
758  flux2 = mh_dh->side_flux( edge_side.side() );
759  if ( flux2 > 0) edg_flux+= flux2;
760  }
761  for( DHCellSide edge_side : cell_side.edge_sides() )
762  if (edge_side != cell_side) {
763  j = edge_side.element().idx();
764  new_j = row_4_el[j];
765 
766  flux2 = mh_dh->side_flux(edge_side.side());
767  if ( flux2 > 0.0 && flux <0.0)
768  aij = -(flux * flux2 / ( edg_flux * dh_cell.elm().measure() ) );
769  else aij =0;
770  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
771  }
772  }
773  if (flux > 0.0)
774  aii -= (flux / dh_cell.elm().measure() );
775  } // end same dim //ELEMENT_SIDES
776 
777  for( DHCellSide neighb_side : dh_cell.neighb_sides() ) // dh_cell lower dim
778  {
779  ASSERT( neighb_side.elem_idx() != dh_cell.elm_idx() ).error("Elm. same\n");
780  new_j = row_4_el[ neighb_side.elem_idx() ];
781  flux = mh_dh->side_flux(neighb_side.side());
782 
783  // volume source - out-flow from higher dimension
784  if (flux > 0.0) aij = flux / dh_cell.elm().measure();
785  else aij=0;
786  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
787  // out flow from higher dim. already accounted
788 
789  // volume drain - in-flow to higher dimension
790  if (flux < 0.0) {
791  aii -= (-flux) / dh_cell.elm().measure(); // diagonal drain
792  aij = (-flux) / neighb_side.element().measure();
793  } else aij=0;
794  MatSetValue(tm, new_j, new_i, aij, INSERT_VALUES);
795  }
796 
797  MatSetValue(tm, new_i, new_i, aii, INSERT_VALUES);
798 
799  cfl_flow_[loc_el++] = fabs(aii);
800  aii = 0.0;
801  }
802 
803  MatAssemblyBegin(tm, MAT_FINAL_ASSEMBLY);
804  MatAssemblyEnd(tm, MAT_FINAL_ASSEMBLY);
805 
807  END_TIMER("convection_matrix_assembly");
808 
810 }
811 
812 
813 
814 
815 
816 //=============================================================================
817 // OUTPUT VECTOR GATHER
818 //=============================================================================
820 
821  unsigned int sbi;
822  IS is;
823 
824  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el, PETSC_COPY_VALUES, &is); //WithArray
825  VecScatterCreate(vconc[0], is, out_conc[0].petsc_vec(), PETSC_NULL, &vconc_out_scatter);
826  for (sbi = 0; sbi < n_substances(); sbi++) {
827  VecScatterBegin(vconc_out_scatter, vconc[sbi], out_conc[sbi].petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
828  VecScatterEnd(vconc_out_scatter, vconc[sbi], out_conc[sbi].petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
829  }
830  chkerr(VecScatterDestroy(&(vconc_out_scatter)));
831  chkerr(ISDestroy(&(is)));
832 }
833 
834 
836  return conc;
837 }
838 
839 void ConvectionTransport::get_par_info(LongIdx * &el_4_loc_out, Distribution * &el_distribution_out){
840  el_4_loc_out = this->el_4_loc;
841  el_distribution_out = this->el_ds;
842  return;
843 }
844 
845 //int *ConvectionTransport::get_el_4_loc(){
846 // return el_4_loc;
847 //}
848 
850  return row_4_el;
851 }
852 
853 
854 
856 
858  //if ( data_.output_fields.is_field_output_time(data_.conc_mobile, time().step()) ) {
860  //}
861 
862  for (unsigned int sbi = 0; sbi < n_substances(); sbi++) {
864  }
865 
866  data_.output_fields.output(time().step());
867 
868  START_TIMER("TOS-balance");
869  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
870  balance_->calculate_instant(sbi, vconc[sbi]);
871  balance_->output();
872  END_TIMER("TOS-balance");
873 }
874 
875 void ConvectionTransport::set_balance_object(std::shared_ptr<Balance> balance)
876 {
877  balance_ = balance;
878  subst_idx = balance_->add_quantities(substances_.names());
879 }
TimeGovernor & time()
Definition: equation.hh:148
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:211
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
unsigned int size() const
get global size
FieldSet * eq_data_
Definition: equation.hh:232
static auto subdomain(Mesh &mesh) -> IndexField
bool is_mass_diag_changed
Flag indicates that porosity or cross_section changed during last time.
Definition: transport.h:268
double time_changed() const
void update_solution() override
Definition: transport.cc:575
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
double end_time() const
End time.
unsigned int * boundary_idx_
Definition: elements.h:83
LongIdx * get_row_4_el() const
Definition: mesh.h:167
double transport_matrix_time
Definition: transport.h:293
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
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
void alloc_transport_vectors()
Definition: transport.cc:274
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:202
void output(TimeStep step)
Boundary * cond() const
Definition: side_impl.hh:71
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:103
void create_mass_matrix()
Definition: transport.cc:700
void next_time()
Proceed to the next time according to current estimated time step.
void initialize() override
Definition: transport.cc:138
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:317
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:275
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:257
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: transport.h:101
double transport_bc_time
Time of the last update of the boundary condition terms.
Definition: transport.h:294
double * cfl_flow_
Definition: transport.h:280
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:361
Definition: mesh.h:76
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:107
Cell accessor allow iterate over DOF handler cells.
double ** sources_corr
Definition: transport.h:271
LongIdx * el_4_loc
Definition: transport.h:321
MultiField< 3, FieldValue< 3 >::Scalar > sources_conc
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_conc
Definition: transport.h:96
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:147
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > output_field_ptr
Fields correspond with out_conc.
Definition: transport.h:312
Definition: edges.h:26
SideIter side(const unsigned int loc_index)
Definition: accessors.hh:137
void fill_output_data(VectorMPI &vec_seq, std::shared_ptr< FieldFE< spacedim, Value > > field_ptr)
Definition: field_fe.hh:309
void set_balance_object(std::shared_ptr< Balance > balance) override
Definition: transport.cc:875
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
double ** conc
Concentrations for phase, substance, element.
Definition: transport.h:299
const TimeStep & step(int index=-1) const
virtual ~ConvectionTransport()
Definition: transport.cc:207
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:99
double t() const
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
void zero_time_step() override
Definition: transport.cc:485
static TimeMarks & marks()
void output_vector_gather()
Definition: transport.cc:819
Basic time management class.
Class for declaration of inputs sequences.
Definition: type_base.hh:346
virtual ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
Definition: mesh.cc:727
const Input::Record input_rec
Record with input specification.
Definition: transport.h:315
static constexpr bool value
Definition: json.hpp:87
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: transport.h:102
const string _equation_name
Definition: transport.cc:61
double last_t() const
const MH_DofHandler * mh_dh
Definition: transport.h:332
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Definition: accessors.hh:285
static const Input::Type::Record & get_input_type()
Definition: transport.cc:67
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void setup_components()
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
RangeConvert< DHEdgeSide, DHCellSide > edge_sides() const
Returns range of all sides looped over common Edge.
FieldCommon & input_default(const string &input_default)
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
unsigned int n_sides() const
Definition: elements.h:135
static auto region_id(Mesh &mesh) -> IndexField
ConvectionTransport(Mesh &init_mesh, const Input::Record in_rec)
Definition: transport.cc:111
static const int registrar
Registrar of class to factory.
Definition: transport.h:252
bool is_changed_dt() const
#define START_TIMER(tag)
Starts a timer with specified tag.
std::vector< VectorMPI > out_conc
Definition: transport.h:307
SubstanceList substances_
Transported substances.
Definition: transport.h:325
Mesh * mesh_
Definition: equation.hh:223
LongIdx * get_row_4_el() override
Return global array of order of elements within parallel vector.
Definition: transport.cc:849
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:389
double ** get_concentration_matrix() override
Getter for array of concentrations per element.
Definition: transport.cc:835
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:503
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1041
std::shared_ptr< DOFHandlerMultiDim > dh_
Definition: transport.h:333
double measure() const
Computes the measure of the element.
Definition: accessors.hh:254
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
Definition: transport.h:336
Distribution * get_el_ds() const
Definition: mesh.h:164
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:297
Region region() const
Definition: accessors.hh:95
void mark_input_times(const TimeGovernor &tg)
Definition: field_set.hh:218
TimeMark add(const TimeMark &mark)
Definition: time_marks.cc:81
void create_transport_matrix_mpi()
Definition: transport.cc:733
Support classes for parallel programing.
void alloc_transport_structs_mpi()
Definition: transport.cc:308
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:265
FLOW123D_FORCE_LINK_IN_CHILD(convectionTransport)
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:190
double * cfl_source_
Definition: transport.h:280
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:346
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
void set_components(const std::vector< string > &names)
Definition: field_set.hh:177
void compute_concentration_sources()
Assembles concentration sources for each substance. note: the source of concentration is multiplied b...
Definition: transport.cc:427
bool evaluate_time_constraint(double &time_constraint) override
Definition: transport.cc:510
VecScatter vconc_out_scatter
Definition: transport.h:283
double dt() const
Vec vcfl_flow_
Parallel vector for flow contribution to CFL condition.
Definition: transport.h:278
void set_target_time(double target_time) override
Definition: transport.cc:677
double ** tm_diag
Definition: transport.h:288
virtual void output_data() override
Write computed fields.
Definition: transport.cc:855
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
void get_par_info(LongIdx *&el_4_loc, Distribution *&el_ds) override
Return array of indices of local elements and parallel distribution of elements.
Definition: transport.cc:839
Distributed sparse graphs, partitioning.
Vec vcfl_source_
Parallel vector for source term contribution to CFL condition.
Definition: transport.h:278
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:91
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
FieldCommon & name(const string &name)
#define END_TIMER(tag)
Ends a timer with specified tag.
#define OLD_ASSERT_EQUAL(a, b)
Definition: global_defs.h:133
double ** cumulative_corr
Definition: transport.h:305
bool changed() const
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:183
Class used for marking specified times at which some events occur.
Definition: time_marks.hh:45
Record type proxy class.
Definition: type_record.hh:182
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:252
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:111
Other possible transformation of coordinates:
LongIdx * get_el_4_loc() const
Definition: mesh.h:170
Distribution * el_ds
Definition: transport.h:322
Implementation of range helper class.
Side accessor allows to iterate over sides of DOF handler cell.
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:48
void make_transport_partitioning()
Definition: transport.cc:180
LongIdx * row_4_el
Definition: transport.h:320
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:276
TimeGovernor * time_
Definition: equation.hh:224
unsigned int lsize(int proc) const
get local size