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