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