Flow123d  JS_before_hm-937-g93502c2
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 #include "system/index_types.hh"
24 
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 "tools/mixed.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 
54 #include "fem/fe_p.hh"
55 #include "fem/fe_values.hh"
57 
58 
59 FLOW123D_FORCE_LINK_IN_CHILD(convectionTransport)
60 
61 
62 namespace IT = Input::Type;
63 
64 /********************************************************************************
65  * Static methods and data members
66  */
67 const string _equation_name = "Solute_Advection_FV";
68 
69 const int ConvectionTransport::registrar =
71  ConvectionTransport::get_input_type().size();
72 
73 const IT::Record &ConvectionTransport::get_input_type()
74 {
75  return IT::Record(_equation_name, "Finite volume method, explicit in time, for advection only solute transport.")
77  .declare_key("input_fields", IT::Array(
78  EqData().make_field_descriptor_type(_equation_name)),
80  "")
81  .declare_key("output",
82  EqData().output_fields.make_output_type(_equation_name, ""),
83  IT::Default("{ \"fields\": [ \"conc\" ] }"),
84  "Specification of output fields and output times.")
85  .close();
86 }
87 
88 
90 {
91  *this += bc_conc.name("bc_conc")
92  .description("Boundary condition for concentration of substances.")
93  .input_default("0.0")
94  .units( UnitSI().kg().m(-3) );
95 
96  *this += init_conc.name("init_conc")
97  .description("Initial values for concentration of substances.")
98  .input_default("0.0")
99  .units( UnitSI().kg().m(-3) );
100 
101  output_fields += *this;
102  output_fields += conc_mobile.name("conc")
103  .units( UnitSI().kg().m(-3) )
105  .description("Concentration solution.");
106  output_fields += region_id.name("region_id")
109  .description("Region ids.");
110  output_fields += subdomain.name("subdomain")
113  .description("Subdomain ids of the domain decomposition.");
114 }
115 
116 
117 /********************************************************************************
118  * Helper class FETransportObjects
119  */
121 : q_(QGauss::make_array(0))
122 {
123  fe0_ = new FE_P_disc<0>(0);
124  fe1_ = new FE_P_disc<1>(0);
125  fe2_ = new FE_P_disc<2>(0);
126  fe3_ = new FE_P_disc<3>(0);
127 
128 
129  fe_values_[0].initialize(q(0), *fe1_,
131  fe_values_[1].initialize(q(1), *fe2_,
133  fe_values_[2].initialize(q(2), *fe3_,
135 }
136 
137 
139 {
140  delete fe0_;
141  delete fe1_;
142  delete fe2_;
143  delete fe3_;
144 }
145 
146 template<> FiniteElement<0> *FETransportObjects::fe<0>() { return fe0_; }
147 template<> FiniteElement<1> *FETransportObjects::fe<1>() { return fe1_; }
148 template<> FiniteElement<2> *FETransportObjects::fe<2>() { return fe2_; }
149 template<> FiniteElement<3> *FETransportObjects::fe<3>() { return fe3_; }
150 
151 Quadrature &FETransportObjects::q(unsigned int dim) { return q_[dim]; }
152 
153 
154 
155 /********************************************************************************
156  * ConvectionTransport
157  */
159 : ConcentrationTransportBase(init_mesh, in_rec),
160  is_mass_diag_changed(false),
161  sources_corr(nullptr),
162  input_rec(in_rec)
163 {
164  START_TIMER("ConvectionTransport");
165  this->eq_data_ = &data_;
166 
167  transport_matrix_time = -1.0; // or -infty
168  transport_bc_time = -1.0;
170  is_src_term_scaled = false;
171  is_bc_term_scaled = false;
172 
173  //initialization of DOF handler
174  MixedPtr<FE_P_disc> fe(0);
175  dh_ = make_shared<DOFHandlerMultiDim>(init_mesh);
176  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>( &init_mesh, fe);
177  dh_->distribute_dofs(ds);
178 
179 }
180 
182 {
184 
186 
188  data_.set_input_list( input_rec.val<Input::Array>("input_fields"), *time_ );
189  data_.set_mesh(*mesh_);
190 
194 
195  // register output vectors
202  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
203  {
204  // create shared pointer to a FieldFE and push this Field to output_field on all regions
205  auto output_field_ptr = create_field_fe< 3, FieldValue<3>::Scalar >(dh_);
206  data_.conc_mobile[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
207  vconc[sbi] = output_field_ptr->get_data_vec().petsc_vec();
208  conc[sbi] = &(output_field_ptr->get_data_vec().data()[0]);
209  }
210  //output_stream_->add_admissible_field_names(input_rec.val<Input::Array>("output_fields"));
211  //output_stream_->mark_output_times(*time_);
212 
214  //cout << "Transport." << endl;
215  //cout << time().marks();
216 
217  balance_->allocate(el_ds->lsize(), 1);
218 
219 }
220 
221 
222 //=============================================================================
223 // MAKE TRANSPORT
224 //=============================================================================
226 
227 // int * id_4_old = new int[mesh_->n_elements()];
228 // int i = 0;
229 // for (auto ele : mesh_->elements_range()) id_4_old[i++] = ele.index();
230 // mesh_->get_part()->id_maps(mesh_->n_elements(), id_4_old, el_ds, el_4_loc, row_4_el);
231 // delete[] id_4_old;
232  el_ds = mesh_->get_el_ds();
235 
236  // TODO: make output of partitioning is usefull but makes outputs different
237  // on different number of processors, which breaks tests.
238  //
239  // Possible solution:
240  // - have flag in ini file to turn this output ON
241  // - possibility to have different ref_output for different num of proc.
242  // - or do not test such kind of output
243  //
244  //for (auto ele : mesh_->elements_range()) {
245  // ele->pid()=el_ds->get_proc(row_4_el[ele.index()]);
246  //}
247 
248 }
249 
250 
251 
253 {
254  unsigned int sbi;
255 
256  if (sources_corr) {
257  //Destroy mpi vectors at first
258  chkerr(MatDestroy(&tm));
259  chkerr(VecDestroy(&mass_diag));
260  chkerr(VecDestroy(&vpmass_diag));
261  chkerr(VecDestroy(&vcfl_flow_));
262  chkerr(VecDestroy(&vcfl_source_));
263  delete cfl_flow_;
264  delete cfl_source_;
265 
266  for (sbi = 0; sbi < n_substances(); sbi++) {
267  // mpi vectors
268  chkerr(VecDestroy(&vpconc[sbi]));
269  chkerr(VecDestroy(&bcvcorr[sbi]));
270  chkerr(VecDestroy(&vcumulative_corr[sbi]));
271  chkerr(VecDestroy(&v_tm_diag[sbi]));
272  chkerr(VecDestroy(&v_sources_corr[sbi]));
273 
274  // arrays of arrays
275  delete cumulative_corr[sbi];
276  delete tm_diag[sbi];
277  delete sources_corr[sbi];
278  }
279 
280  // arrays of mpi vectors
281  delete vpconc;
282  delete bcvcorr;
283  delete vcumulative_corr;
284  delete v_tm_diag;
285  delete v_sources_corr;
286 
287  // arrays of arrays
288  delete conc;
289  delete cumulative_corr;
290  delete tm_diag;
291  delete sources_corr;
292  }
293 }
294 
295 
296 
297 
298 
300 {
301  for ( DHCellAccessor dh_cell : dh_->own_range() ) {
302  LongIdx index = dh_cell.local_idx();
303  ElementAccessor<3> ele_acc = mesh_->element_accessor( dh_cell.elm_idx() );
304 
305  for (unsigned int sbi=0; sbi<n_substances(); sbi++) // Optimize: SWAP LOOPS
306  conc[sbi][index] = data_.init_conc[sbi].value(ele_acc.centre(), ele_acc);
307  }
308 }
309 
310 //=============================================================================
311 // ALLOCATE OF TRANSPORT VARIABLES (ELEMENT & NODES)
312 //=============================================================================
314 
315  unsigned int sbi, n_subst;
316  n_subst = n_substances();
317 
318  sources_corr = new double*[n_subst];
319  tm_diag = new double*[n_subst];
320  cumulative_corr = new double*[n_subst];
321  for (sbi = 0; sbi < n_subst; sbi++) {
322  cumulative_corr[sbi] = new double[el_ds->lsize()];
323  sources_corr[sbi] = new double[el_ds->lsize()];
324  tm_diag[sbi] = new double[el_ds->lsize()];
325  }
326 
327  conc = new double*[n_subst];
328  vconc = new Vec[n_subst];
329 
330  cfl_flow_ = new double[el_ds->lsize()];
331  cfl_source_ = new double[el_ds->lsize()];
332 }
333 
334 //=============================================================================
335 // ALLOCATION OF TRANSPORT VECTORS (MPI)
336 //=============================================================================
338 
339  int sbi, n_subst;
340  n_subst = n_substances();
341 
342  MPI_Barrier(PETSC_COMM_WORLD);
343 
344  vpconc = new Vec[n_subst];
345  bcvcorr = new Vec[n_subst];
346  vcumulative_corr = new Vec[n_subst];
347  v_tm_diag = new Vec[n_subst];
348  v_sources_corr = new Vec[n_subst];
349 
350 
351  for (sbi = 0; sbi < n_subst; sbi++) {
352  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &bcvcorr[sbi]);
353  VecZeroEntries(bcvcorr[sbi]);
354 
355  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpconc[sbi]);
356  VecZeroEntries(vpconc[sbi]);
357 
358  // SOURCES
359  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
360  cumulative_corr[sbi],&vcumulative_corr[sbi]);
361 
362  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
363  sources_corr[sbi],&v_sources_corr[sbi]);
364 
365  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
366  tm_diag[sbi],&v_tm_diag[sbi]);
367 
368  VecZeroEntries(vcumulative_corr[sbi]);
369  }
370 
371 
372  MatCreateAIJ(PETSC_COMM_WORLD, el_ds->lsize(), el_ds->lsize(), mesh_->n_elements(),
373  mesh_->n_elements(), 16, PETSC_NULL, 4, PETSC_NULL, &tm);
374 
375  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &mass_diag);
376  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpmass_diag);
377 
378  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
380  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
382 }
383 
384 
386 {
387  START_TIMER ("set_boundary_conditions");
388 
389  unsigned int sbi;
390 
391  // Assembly bcvcorr vector
392  for(sbi=0; sbi < n_substances(); sbi++) VecZeroEntries(bcvcorr[sbi]);
393 
394  balance_->start_flux_assembly(subst_idx);
395 
396  for ( DHCellAccessor dh_cell : dh_->own_range() )
397  {
398  ElementAccessor<3> elm = dh_cell.elm();
399  // we have currently zero order P_Disc FE
400  ASSERT_DBG(dh_cell.get_loc_dof_indices().size() == 1);
401  IntIdx local_p0_dof = dh_cell.get_loc_dof_indices()[0];
402  LongIdx glob_p0_dof = dh_->get_local_to_global_map()[local_p0_dof];
403 
404  for(DHCellSide dh_side: dh_cell.side_range()) {
405  if (dh_side.side().is_boundary()) {
406  ElementAccessor<3> bc_elm = dh_side.cond().element_accessor();
407  double flux = this->side_flux(dh_side);
408  if (flux < 0.0) {
409  double aij = -(flux / elm.measure() );
410 
411  for (sbi=0; sbi<n_substances(); sbi++)
412  {
413  double value = data_.bc_conc[sbi].value( bc_elm.centre(), bc_elm );
414 
415  VecSetValue(bcvcorr[sbi], glob_p0_dof, value * aij, ADD_VALUES);
416 
417  // CAUTION: It seems that PETSc possibly optimize allocated space during assembly.
418  // So we have to add also values that may be non-zero in future due to changing velocity field.
419  balance_->add_flux_values(subst_idx[sbi], dh_side,
420  {local_p0_dof}, {0.0}, flux*value);
421  }
422  } else {
423  for (sbi=0; sbi<n_substances(); sbi++)
424  VecSetValue(bcvcorr[sbi], glob_p0_dof, 0, ADD_VALUES);
425 
426  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
427  balance_->add_flux_values(subst_idx[sbi], dh_side,
428  {local_p0_dof}, {flux}, 0.0);
429  }
430  }
431  }
432  }
433 
434  balance_->finish_flux_assembly(subst_idx);
435 
436  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyBegin(bcvcorr[sbi]);
437  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyEnd(bcvcorr[sbi]);
438 
439  // we are calling set_boundary_conditions() after next_time() and
440  // we are using data from t() before, so we need to set corresponding bc time
442 }
443 
444 
445 //=============================================================================
446 // COMPUTE SOURCES
447 //=============================================================================
449 
450  //temporary variables
451  double csection, source, diag;
452 
453  //TODO: would it be possible to check the change in data for chosen substance? (may be in multifields?)
454 
455  //checking if the data were changed
456  if( (data_.sources_density.changed() )
457  || (data_.sources_conc.changed() )
458  || (data_.sources_sigma.changed() )
459  || (data_.cross_section.changed()))
460  {
461  START_TIMER("sources_reinit");
462  balance_->start_source_assembly(subst_idx);
463 
464  for ( DHCellAccessor dh_cell : dh_->own_range() )
465  {
466  ElementAccessor<3> elm = dh_cell.elm();
467  // we have currently zero order P_Disc FE
468  ASSERT_DBG(dh_cell.get_loc_dof_indices().size() == 1);
469  IntIdx local_p0_dof = dh_cell.get_loc_dof_indices()[0];
470 
471  arma::vec3 center = elm.centre();
472  csection = data_.cross_section.value(center, elm);
473 
474  // read for all substances
475  double max_cfl=0;
476  for (unsigned int sbi = 0; sbi < n_substances(); sbi++)
477  {
478  double src_sigma = data_.sources_sigma[sbi].value(center, elm);
479 
480  source = csection * (data_.sources_density[sbi].value(center, elm)
481  + src_sigma * data_.sources_conc[sbi].value(center, elm));
482  // addition to RHS
483  sources_corr[sbi][local_p0_dof] = source;
484  // addition to diagonal of the transport matrix
485  diag = src_sigma * csection;
486  tm_diag[sbi][local_p0_dof] = - diag;
487 
488  // compute maximal cfl condition over all substances
489  max_cfl = std::max(max_cfl, fabs(diag));
490 
491  balance_->add_source_values(sbi, elm.region().bulk_idx(), {local_p0_dof},
492  {- src_sigma * elm.measure() * csection},
493  {source * elm.measure()});
494  }
495 
496  cfl_source_[local_p0_dof] = max_cfl;
497  }
498 
499  balance_->finish_source_assembly(subst_idx);
500 
501  END_TIMER("sources_reinit");
502  }
503 }
504 
505 
506 
508 {
510 
513  std::stringstream ss; // print warning message with table of uninitialized fields
514  if ( FieldCommon::print_message_table(ss, "convection transport") ) {
515  WarningOut() << ss.str();
516  }
517 
520 
521  START_TIMER("Convection balance zero time step");
522 
526 
527  // write initial condition
528  output_data();
529 }
530 
531 
533 {
534  ASSERT_PTR(dh_).error( "Null DOF handler object.\n" );
535  // read changed status before setting time
536  bool changed_flux = data_.flow_flux.changed();
537  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
538 
539  START_TIMER("data reinit");
540 
541  bool cfl_changed = false;
542 
543  // if FLOW or DATA changed ---------------------> recompute transport matrix
544  if (changed_flux)
545  {
548  cfl_changed = true;
549  DebugOut() << "CFL changed - flow.\n";
550  }
551 
553  {
555  cfl_changed = true;
556  DebugOut() << "CFL changed - mass matrix.\n";
557  }
558 
559  // if DATA changed ---------------------> recompute concentration sources (rhs and matrix diagonal)
562  {
564  is_src_term_scaled = false;
565  cfl_changed = true;
566  DebugOut() << "CFL changed - source.\n";
567  }
568 
569  // now resolve the CFL condition
570  if(cfl_changed)
571  {
572  // find maximum of sum of contribution from flow and sources: MAX(vcfl_flow_ + vcfl_source_)
573  Vec cfl;
574  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(),PETSC_DETERMINE, &cfl);
575  VecWAXPY(cfl, 1.0, vcfl_flow_, vcfl_source_);
576  VecMaxPointwiseDivide(cfl,mass_diag, &cfl_max_step);
577  // get a reciprocal value as a time constraint
579  DebugOut().fmt("CFL constraint (transport): {}\n", cfl_max_step);
580  }
581 
582  // although it does not influence CFL, compute BC so the full system is assembled
583  if ( data_.flow_flux.changed()
584  || data_.porosity.changed()
586  || data_.bc_conc.changed() )
587  {
589  is_bc_term_scaled = false;
590  }
591 
592  END_TIMER("data reinit");
593 
594  // return time constraint
595  time_constraint = cfl_max_step;
596  return cfl_changed;
597 }
598 
600 
601  START_TIMER("convection-one step");
602 
603  // proceed to next time - which we are about to compute
604  // explicit scheme looks one step back and uses data from previous time
605  // (data time set previously in assess_time_constraint())
606  time_->next_time();
607 
608  double dt_new = time_->dt(), // current time step we are about to compute
609  dt_scaled = dt_new / time_->last_dt(); // scaling ratio to previous time step
610 
611  START_TIMER("time step rescaling");
612 
613  // if FLOW or DATA or BC or DT changed ---------------------> rescale boundary condition
615  {
616  DebugOut() << "BC - rescale dt.\n";
617  //choose between fresh scaling with new dt or rescaling to a new dt
618  double dt = (!is_bc_term_scaled) ? dt_new : dt_scaled;
619  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
620  VecScale(bcvcorr[sbi], dt);
621  is_bc_term_scaled = true;
622  }
623 
624 
625  // if DATA or TIME STEP changed -----------------------> rescale source term
627  DebugOut() << "SRC - rescale dt.\n";
628  //choose between fresh scaling with new dt or rescaling to a new dt
629  double dt = (!is_src_term_scaled) ? dt_new : dt_scaled;
630  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
631  {
632  VecScale(v_sources_corr[sbi], dt);
633  VecScale(v_tm_diag[sbi], dt);
634  }
635  is_src_term_scaled = true;
636  }
637 
638  // if DATA or TIME STEP changed -----------------------> rescale transport matrix
640  DebugOut() << "TM - rescale dt.\n";
641  //choose between fresh scaling with new dt or rescaling to a new dt
642  double dt = (!is_convection_matrix_scaled) ? dt_new : dt_scaled;
643 
644  MatScale(tm, dt);
646  }
647 
648  END_TIMER("time step rescaling");
649 
650 
651  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
653  {
654  VecCopy(mass_diag, vpmass_diag);
656  } else is_mass_diag_changed = false;
657 
658 
659  // Compute new concentrations for every substance.
660 
661  for (unsigned int sbi = 0; sbi < n_substances(); sbi++) {
662  // one step in MOBILE phase
663  START_TIMER("mat mult");
664 
665  // tm_diag is a diagonal part of transport matrix, which depends on substance data (sources_sigma)
666  // Wwe need keep transport matrix independent of substance, therefore we keep this diagonal part
667  // separately in a vector tm_diag.
668  // Computation: first, we compute this diagonal addition D*pconc and save it temporaly into RHS
669 
670  // RHS = D*pconc, where D is diagonal matrix represented by a vector
671  VecPointwiseMult(vcumulative_corr[sbi], v_tm_diag[sbi], vconc[sbi]); //w = x.*y
672 
673  // Then we add boundary terms ans other source terms into RHS.
674  // RHS = 1.0 * bcvcorr + 1.0 * v_sources_corr + 1.0 * rhs
675  VecAXPBYPCZ(vcumulative_corr[sbi], 1.0, 1.0, 1.0, bcvcorr[sbi], v_sources_corr[sbi]); //z = ax + by + cz
676 
677  // Then we set the new previous concentration.
678  VecCopy(vconc[sbi], vpconc[sbi]); // pconc = conc
679  // And finally proceed with transport matrix multiplication.
680  if (is_mass_diag_changed) {
681  VecPointwiseMult(vconc[sbi], vconc[sbi], vpmass_diag); // vconc*=vpmass_diag
682  MatMultAdd(tm, vpconc[sbi], vconc[sbi], vconc[sbi]); // vconc+=tm*vpconc
683  VecAXPY(vconc[sbi], 1, vcumulative_corr[sbi]); // vconc+=vcumulative_corr
684  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
685  } else {
686  MatMultAdd(tm, vpconc[sbi], vcumulative_corr[sbi], vconc[sbi]); // vconc =tm*vpconc+vcumulative_corr
687  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
688  VecAXPY(vconc[sbi], 1, vpconc[sbi]); // vconc+=vpconc
689  }
690 
691  END_TIMER("mat mult");
692  }
693 
694  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
695  balance_->calculate_cumulative(sbi, vpconc[sbi]);
696 
697  END_TIMER("convection-one step");
698 }
699 
700 
701 void ConvectionTransport::set_target_time(double target_time)
702 {
703 
704  //sets target_mark_type (it is fixed) to be met in next_time()
705  time_->marks().add(TimeMark(target_time, target_mark_type));
706 
707  // This is done every time TOS calls update_solution.
708  // If CFL condition is changed, time fixation will change later from TOS.
709 
710  // Set the same constraint as was set last time.
711 
712  // TODO: fix this hack, remove this method completely, leaving just the first line at the calling point
713  // in TransportOperatorSplitting::update_solution()
714  // doing this directly leads to choose of large time step violationg CFL condition
715  if (cfl_max_step > time_->dt()*1e-10)
716  time_->set_upper_constraint(cfl_max_step, "CFL condition used from previous step.");
717 
718  // fixing convection time governor till next target_mark_type (got from TOS or other)
719  // may have marks for data changes
721 }
722 
723 
725 {
726  ElementAccessor<3> elm;
727 
728  VecZeroEntries(mass_diag);
729 
730  balance_->start_mass_assembly(subst_idx);
731 
732  for ( DHCellAccessor dh_cell : dh_->own_range() ) {
733  ElementAccessor<3> elm = dh_cell.elm();
734  // we have currently zero order P_Disc FE
735  ASSERT_DBG(dh_cell.get_loc_dof_indices().size() == 1);
736  IntIdx local_p0_dof = dh_cell.get_loc_dof_indices()[0];
737 
738  double csection = data_.cross_section.value(elm.centre(), elm);
739  //double por_m = data_.porosity.value(elm.centre(), elm->element_accessor());
740  double por_m = data_.water_content.value(elm.centre(), elm);
741 
742  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
743  balance_->add_mass_values(subst_idx[sbi], dh_cell, {local_p0_dof}, {csection*por_m*elm.measure()}, 0);
744 
745  VecSetValue(mass_diag, dh_->get_local_to_global_map()[local_p0_dof], csection*por_m, INSERT_VALUES);
746  }
747 
748  balance_->finish_mass_assembly(subst_idx);
749 
750  VecAssemblyBegin(mass_diag);
751  VecAssemblyEnd(mass_diag);
752 
753  is_mass_diag_changed = true;
754 }
755 
756 
757 //=============================================================================
758 // CREATE TRANSPORT MATRIX
759 //=============================================================================
761 
762  START_TIMER("convection_matrix_assembly");
763 
764  ElementAccessor<3> el2;
765  ElementAccessor<3> elm;
766  int j;
767  LongIdx new_j, new_i;
768  double aij, aii;
769 
770  MatZeroEntries(tm);
771 
772  double flux, flux2, edg_flux;
773 
774  aii = 0.0;
775 
776  unsigned int loc_el = 0;
777  for ( DHCellAccessor dh_cell : dh_->own_range() ) {
778  new_i = row_4_el[ dh_cell.elm_idx() ];
779  elm = dh_cell.elm();
780  for( DHCellSide cell_side : dh_cell.side_range() ) {
781  flux = this->side_flux(cell_side);
782  if (! cell_side.side().is_boundary()) {
783  edg_flux = 0;
784  for( DHCellSide edge_side : cell_side.edge_sides() ) {
785  el2 = edge_side.element();
786  flux2 = this->side_flux(edge_side);
787  if ( flux2 > 0) edg_flux+= flux2;
788  }
789  for( DHCellSide edge_side : cell_side.edge_sides() )
790  if (edge_side != cell_side) {
791  j = edge_side.element().idx();
792  new_j = row_4_el[j];
793 
794  el2 = edge_side.element();
795  flux2 = this->side_flux(edge_side);
796  if ( flux2 > 0.0 && flux <0.0)
797  aij = -(flux * flux2 / ( edg_flux * dh_cell.elm().measure() ) );
798  else aij =0;
799  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
800  }
801  }
802  if (flux > 0.0)
803  aii -= (flux / dh_cell.elm().measure() );
804  } // end same dim //ELEMENT_SIDES
805 
806  for( DHCellSide neighb_side : dh_cell.neighb_sides() ) // dh_cell lower dim
807  {
808  ASSERT( neighb_side.elem_idx() != dh_cell.elm_idx() ).error("Elm. same\n");
809  new_j = row_4_el[ neighb_side.elem_idx() ];
810  el2 = neighb_side.element();
811  flux = this->side_flux(neighb_side);
812 
813  // volume source - out-flow from higher dimension
814  if (flux > 0.0) aij = flux / dh_cell.elm().measure();
815  else aij=0;
816  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
817  // out flow from higher dim. already accounted
818 
819  // volume drain - in-flow to higher dimension
820  if (flux < 0.0) {
821  aii -= (-flux) / dh_cell.elm().measure(); // diagonal drain
822  aij = (-flux) / neighb_side.element().measure();
823  } else aij=0;
824  MatSetValue(tm, new_j, new_i, aij, INSERT_VALUES);
825  }
826 
827  MatSetValue(tm, new_i, new_i, aii, INSERT_VALUES);
828 
829  cfl_flow_[loc_el++] = fabs(aii);
830  aii = 0.0;
831  }
832 
833  MatAssemblyBegin(tm, MAT_FINAL_ASSEMBLY);
834  MatAssemblyEnd(tm, MAT_FINAL_ASSEMBLY);
835 
837  END_TIMER("convection_matrix_assembly");
838 
840 }
841 
842 
844  return conc;
845 }
846 
847 void ConvectionTransport::get_par_info(LongIdx * &el_4_loc_out, Distribution * &el_distribution_out){
848  el_4_loc_out = this->el_4_loc;
849  el_distribution_out = this->el_ds;
850  return;
851 }
852 
853 //int *ConvectionTransport::get_el_4_loc(){
854 // return el_4_loc;
855 //}
856 
858  return row_4_el;
859 }
860 
861 
862 
864 
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 }
880 
882 {
883  if (cell_side.dim()==3) return calculate_side_flux<3>(cell_side);
884  else if (cell_side.dim()==2) return calculate_side_flux<2>(cell_side);
885  else return calculate_side_flux<1>(cell_side);
886 }
887 
888 template<unsigned int dim>
890 {
891  ASSERT_EQ(cell_side.dim(), dim).error("Element dimension mismatch!");
892 
893  feo_.fe_values(dim).reinit(cell_side.side());
894  auto vel = data_.flow_flux.value(cell_side.centre(), cell_side.element());
895  double side_flux = arma::dot(vel, feo_.fe_values(dim).normal_vector(0)) * feo_.fe_values(dim).JxW(0);
896  return side_flux;
897 }
TimeGovernor & time()
Definition: equation.hh:151
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:214
Shape function values.
Definition: update_flags.hh:87
FieldSet * eq_data_
Definition: equation.hh:229
const Element * element() const
Definition: accessors.hh:160
static auto subdomain(Mesh &mesh) -> IndexField
Quadrature & q(unsigned int dim)
Definition: transport.cc:151
bool is_mass_diag_changed
Flag indicates that porosity or cross_section changed during last time.
Definition: transport.h:310
Transformed quadrature weight for cell sides.
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:329
void update_solution() override
Definition: transport.cc:599
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
double end_time() const
End time.
LongIdx * get_row_4_el() const
Definition: mesh.h:157
double transport_matrix_time
Definition: transport.h:335
Field< 3, FieldValue< 3 >::VectorFixed > flow_flux
Flow flux, can be result of water flow model.
int tlevel() const
double calculate_side_flux(const DHCellSide &cell)
Calculate flux on side of given element specified by dimension.
Definition: transport.cc:889
unsigned int dim() const
Return dimension of element appropriate to the side.
const std::vector< std::string > & names()
Definition: substance.hh:85
QGauss::array q_
Quadratures used in assembling methods.
Definition: transport.h:107
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:313
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:238
int IntIdx
Definition: index_types.hh:25
void output(TimeStep step)
double fix_dt_until_mark()
Fixing time step until fixed time mark.
Abstract linear system class.
Definition: balance.hh:40
MultiField< 3, FieldValue< 3 >::Scalar > conc_mobile
Calculated concentrations in the mobile zone.
Definition: transport.h:141
void create_mass_matrix()
Definition: transport.cc:724
void next_time()
Proceed to the next time according to current estimated time step.
void initialize() override
Definition: transport.cc:181
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:352
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:317
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:299
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: transport.h:139
double transport_bc_time
Time of the last update of the boundary condition terms.
Definition: transport.h:336
double * cfl_flow_
Definition: transport.h:322
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:385
Definition: mesh.h:78
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:145
Cell accessor allow iterate over DOF handler cells.
double ** sources_corr
Definition: transport.h:313
LongIdx * el_4_loc
Definition: transport.h:356
Class FEValues calculates finite element data on the actual cells such as shape function values...
MultiField< 3, FieldValue< 3 >::Scalar > sources_conc
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_conc
Definition: transport.h:134
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:148
double side_flux(const DHCellSide &cell_side)
Wrapper of following method, call side_flux with correct template parameter.
Definition: transport.cc:881
const RegionDB & region_db() const
Definition: mesh.h:143
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:347
double ** conc
Concentrations for phase, substance, element.
Definition: transport.h:341
const TimeStep & step(int index=-1) const
virtual ~ConvectionTransport()
Definition: transport.cc:252
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:137
double t() const
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
void zero_time_step() override
Definition: transport.cc:507
Side side() const
Return Side of given cell and side_idx.
int register_class(string class_name)
Function allows simplified call of registering class to factory.
Definition: factory_impl.hh:64
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:548
static TimeMarks & marks()
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:839
const Input::Record input_rec
Record with input specification.
Definition: transport.h:350
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
static constexpr bool value
Definition: json.hpp:87
Symmetric Gauss-Legendre quadrature formulae on simplices.
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: transport.h:140
const string _equation_name
Definition: transport.cc:67
double last_t() const
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
void setup_components()
Discontinuous Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:110
Transformed quadrature points.
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) ...
Definitions of basic Lagrangean finite elements with polynomial shape functions.
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:58
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
ElementAccessor< 3 > element() const
const Ret val(const string &key) const
static auto region_id(Mesh &mesh) -> IndexField
ConvectionTransport(Mesh &init_mesh, const Input::Record in_rec)
Definition: transport.cc:158
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
Definition: fe_values.cc:136
bool is_changed_dt() const
#define START_TIMER(tag)
Starts a timer with specified tag.
SubstanceList substances_
Transported substances.
Definition: transport.h:360
Mesh * mesh_
Definition: equation.hh:220
LongIdx * get_row_4_el() override
Return global array of order of elements within parallel vector.
Definition: transport.cc:857
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:434
FEValues< 3 > & fe_values(unsigned int dim)
Definition: transport.h:92
double ** get_concentration_matrix() override
Getter for array of concentrations per element.
Definition: transport.cc:843
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:232
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
std::shared_ptr< DOFHandlerMultiDim > dh_
Definition: transport.h:367
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:370
Shape function gradients.
Definition: update_flags.hh:95
Distribution * get_el_ds() const
Definition: mesh.h:154
std::shared_ptr< Balance > balance() const
Definition: equation.hh:187
FEValues< 3 > fe_values_[3]
FESideValues objects for side flux calculating.
Definition: transport.h:110
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:339
Region region() const
Definition: accessors.hh:165
Normal vectors.
void mark_input_times(const TimeGovernor &tg)
Definition: field_set.hh:221
TimeMark add(const TimeMark &mark)
Definition: time_marks.cc:81
void create_transport_matrix_mpi()
Definition: transport.cc:760
Support classes for parallel programing.
void alloc_transport_structs_mpi()
Definition: transport.cc:337
int set_upper_constraint(double upper, std::string message)
Sets upper constraint for the next time step estimating.
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
FieldCommon & description(const string &description)
bool is_convection_matrix_scaled
Definition: transport.h:307
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:193
double * cfl_source_
Definition: transport.h:322
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:337
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
void set_components(const std::vector< string > &names)
Definition: field_set.hh:180
void compute_concentration_sources()
Assembles concentration sources for each substance. note: the source of concentration is multiplied b...
Definition: transport.cc:448
bool evaluate_time_constraint(double &time_constraint) override
Definition: transport.cc:532
FiniteElement< 1 > * fe1_
Definition: transport.h:102
double dt() const
Vec vcfl_flow_
Parallel vector for flow contribution to CFL condition.
Definition: transport.h:320
void set_target_time(double target_time) override
Definition: transport.cc:701
double ** tm_diag
Definition: transport.h:330
virtual void output_data() override
Write computed fields.
Definition: transport.cc:863
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:847
Distributed sparse graphs, partitioning.
#define ASSERT_DBG(expr)
Vec vcfl_source_
Parallel vector for source term contribution to CFL condition.
Definition: transport.h:320
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:91
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:258
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:347
bool changed() const
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:186
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
FETransportObjects feo_
Finite element objects.
Definition: transport.h:373
FieldCommon & flags(FieldFlag::Flags::Mask mask)
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.hh:200
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:264
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
Other possible transformation of coordinates:
FiniteElement< 3 > * fe3_
Definition: transport.h:104
FiniteElement< 0 > * fe0_
Finite elements for the solution of the advection-diffusion equation.
Definition: transport.h:101
LongIdx * get_el_4_loc() const
Definition: mesh.h:160
Distribution * el_ds
Definition: transport.h:357
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
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
void make_transport_partitioning()
Definition: transport.cc:225
LongIdx * row_4_el
Definition: transport.h:355
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:318
TimeGovernor * time_
Definition: equation.hh:221
FiniteElement< 2 > * fe2_
Definition: transport.h:103
unsigned int lsize(int proc) const
get local size
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:234
arma::vec3 centre() const
Side centre.