Flow123d  JS_before_hm-929-gaeebe69
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  changed_(true)
164 {
165  START_TIMER("ConvectionTransport");
166  this->eq_data_ = &data_;
167 
168  transport_matrix_time = -1.0; // or -infty
169  transport_bc_time = -1.0;
171  is_src_term_scaled = false;
172  is_bc_term_scaled = false;
173 
174  //initialization of DOF handler
175  MixedPtr<FE_P_disc> fe(0);
176  dh_ = make_shared<DOFHandlerMultiDim>(init_mesh);
177  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>( &init_mesh, fe);
178  dh_->distribute_dofs(ds);
179 
180 }
181 
183 {
185 
187 
189  data_.set_input_list( input_rec.val<Input::Array>("input_fields"), *time_ );
190  data_.set_mesh(*mesh_);
191 
195 
196  // register output vectors
203  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
204  {
205  // create shared pointer to a FieldFE and push this Field to output_field on all regions
206  output_field_ptr[sbi] = std::make_shared< FieldFE<3, FieldValue<3>::Scalar> >();
207  output_field_ptr[sbi]->set_fe_data(dh_ );
208  data_.conc_mobile[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr[sbi], 0);
209  vconc[sbi] = output_field_ptr[sbi]->get_data_vec().petsc_vec();
210  conc[sbi] = &(output_field_ptr[sbi]->get_data_vec().data()[0]);
211  }
212  //output_stream_->add_admissible_field_names(input_rec.val<Input::Array>("output_fields"));
213  //output_stream_->mark_output_times(*time_);
214 
216  //cout << "Transport." << endl;
217  //cout << time().marks();
218 
219  balance_->allocate(el_ds->lsize(), 1);
220 
221 }
222 
223 
224 //=============================================================================
225 // MAKE TRANSPORT
226 //=============================================================================
228 
229 // int * id_4_old = new int[mesh_->n_elements()];
230 // int i = 0;
231 // for (auto ele : mesh_->elements_range()) id_4_old[i++] = ele.index();
232 // mesh_->get_part()->id_maps(mesh_->n_elements(), id_4_old, el_ds, el_4_loc, row_4_el);
233 // delete[] id_4_old;
234  el_ds = mesh_->get_el_ds();
237 
238  // TODO: make output of partitioning is usefull but makes outputs different
239  // on different number of processors, which breaks tests.
240  //
241  // Possible solution:
242  // - have flag in ini file to turn this output ON
243  // - possibility to have different ref_output for different num of proc.
244  // - or do not test such kind of output
245  //
246  //for (auto ele : mesh_->elements_range()) {
247  // ele->pid()=el_ds->get_proc(row_4_el[ele.index()]);
248  //}
249 
250 }
251 
252 
253 
255 {
256  unsigned int sbi;
257 
258  if (sources_corr) {
259  //Destroy mpi vectors at first
260  chkerr(MatDestroy(&tm));
261  chkerr(VecDestroy(&mass_diag));
262  chkerr(VecDestroy(&vpmass_diag));
263  chkerr(VecDestroy(&vcfl_flow_));
264  chkerr(VecDestroy(&vcfl_source_));
265  delete cfl_flow_;
266  delete cfl_source_;
267 
268  for (sbi = 0; sbi < n_substances(); sbi++) {
269  // mpi vectors
270  chkerr(VecDestroy(&vpconc[sbi]));
271  chkerr(VecDestroy(&bcvcorr[sbi]));
272  chkerr(VecDestroy(&vcumulative_corr[sbi]));
273  chkerr(VecDestroy(&v_tm_diag[sbi]));
274  chkerr(VecDestroy(&v_sources_corr[sbi]));
275 
276  // arrays of arrays
277  delete cumulative_corr[sbi];
278  delete tm_diag[sbi];
279  delete sources_corr[sbi];
280  }
281 
282  // arrays of mpi vectors
283  delete vpconc;
284  delete bcvcorr;
285  delete vcumulative_corr;
286  delete v_tm_diag;
287  delete v_sources_corr;
288 
289  // arrays of arrays
290  delete conc;
291  delete cumulative_corr;
292  delete tm_diag;
293  delete sources_corr;
294  }
295 }
296 
297 
298 
299 
300 
302 {
303  for ( DHCellAccessor dh_cell : dh_->own_range() ) {
304  LongIdx index = dh_cell.local_idx();
305  ElementAccessor<3> ele_acc = mesh_->element_accessor( dh_cell.elm_idx() );
306 
307  for (unsigned int sbi=0; sbi<n_substances(); sbi++) // Optimize: SWAP LOOPS
308  output_field_ptr[sbi]->get_data_vec().data()[index] = data_.init_conc[sbi].value(ele_acc.centre(), ele_acc);
309  }
310 }
311 
312 //=============================================================================
313 // ALLOCATE OF TRANSPORT VARIABLES (ELEMENT & NODES)
314 //=============================================================================
316 
317  unsigned int sbi, n_subst;
318  n_subst = n_substances();
319 
320  sources_corr = new double*[n_subst];
321  tm_diag = new double*[n_subst];
322  cumulative_corr = new double*[n_subst];
323  for (sbi = 0; sbi < n_subst; sbi++) {
324  cumulative_corr[sbi] = new double[el_ds->lsize()];
325  sources_corr[sbi] = new double[el_ds->lsize()];
326  tm_diag[sbi] = new double[el_ds->lsize()];
327  }
328 
329  conc = new double*[n_subst];
330  vconc = new Vec[n_subst];
331  out_conc.clear();
332  out_conc.resize(n_subst);
333  output_field_ptr.clear();
334  output_field_ptr.resize(n_subst);
335  for (sbi = 0; sbi < n_subst; sbi++) {
336  out_conc[sbi].resize( el_ds->size() );
337  }
338 
339  cfl_flow_ = new double[el_ds->lsize()];
340  cfl_source_ = new double[el_ds->lsize()];
341 }
342 
343 //=============================================================================
344 // ALLOCATION OF TRANSPORT VECTORS (MPI)
345 //=============================================================================
347 
348  int sbi, n_subst;
349  n_subst = n_substances();
350 
351  MPI_Barrier(PETSC_COMM_WORLD);
352 
353  vpconc = new Vec[n_subst];
354  bcvcorr = new Vec[n_subst];
355  vcumulative_corr = new Vec[n_subst];
356  v_tm_diag = new Vec[n_subst];
357  v_sources_corr = new Vec[n_subst];
358 
359 
360  for (sbi = 0; sbi < n_subst; sbi++) {
361  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &bcvcorr[sbi]);
362  VecZeroEntries(bcvcorr[sbi]);
363 
364  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpconc[sbi]);
365  VecZeroEntries(vpconc[sbi]);
366 
367  // SOURCES
368  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
369  cumulative_corr[sbi],&vcumulative_corr[sbi]);
370 
371  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
372  sources_corr[sbi],&v_sources_corr[sbi]);
373 
374  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
375  tm_diag[sbi],&v_tm_diag[sbi]);
376 
377  VecZeroEntries(vcumulative_corr[sbi]);
378  VecZeroEntries(out_conc[sbi].petsc_vec());
379  }
380 
381 
382  MatCreateAIJ(PETSC_COMM_WORLD, el_ds->lsize(), el_ds->lsize(), mesh_->n_elements(),
383  mesh_->n_elements(), 16, PETSC_NULL, 4, PETSC_NULL, &tm);
384 
385  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &mass_diag);
386  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpmass_diag);
387 
388  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
390  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
392 }
393 
394 
396 {
397  START_TIMER ("set_boundary_conditions");
398 
399  unsigned int sbi;
400 
401  // Assembly bcvcorr vector
402  for(sbi=0; sbi < n_substances(); sbi++) VecZeroEntries(bcvcorr[sbi]);
403 
404  balance_->start_flux_assembly(subst_idx);
405 
406  for ( DHCellAccessor dh_cell : dh_->own_range() )
407  {
408  ElementAccessor<3> elm = dh_cell.elm();
409  // we have currently zero order P_Disc FE
410  ASSERT_DBG(dh_cell.get_loc_dof_indices().size() == 1);
411  IntIdx local_p0_dof = dh_cell.get_loc_dof_indices()[0];
412  LongIdx glob_p0_dof = dh_->get_local_to_global_map()[local_p0_dof];
413 
414  for(DHCellSide dh_side: dh_cell.side_range()) {
415  if (dh_side.side().is_boundary()) {
416  ElementAccessor<3> bc_elm = dh_side.cond().element_accessor();
417  double flux = this->side_flux(dh_side);
418  if (flux < 0.0) {
419  double aij = -(flux / elm.measure() );
420 
421  for (sbi=0; sbi<n_substances(); sbi++)
422  {
423  double value = data_.bc_conc[sbi].value( bc_elm.centre(), bc_elm );
424 
425  VecSetValue(bcvcorr[sbi], glob_p0_dof, value * aij, ADD_VALUES);
426 
427  // CAUTION: It seems that PETSc possibly optimize allocated space during assembly.
428  // So we have to add also values that may be non-zero in future due to changing velocity field.
429  balance_->add_flux_values(subst_idx[sbi], dh_side,
430  {local_p0_dof}, {0.0}, flux*value);
431  }
432  } else {
433  for (sbi=0; sbi<n_substances(); sbi++)
434  VecSetValue(bcvcorr[sbi], glob_p0_dof, 0, ADD_VALUES);
435 
436  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
437  balance_->add_flux_values(subst_idx[sbi], dh_side,
438  {local_p0_dof}, {flux}, 0.0);
439  }
440  }
441  }
442  }
443 
444  balance_->finish_flux_assembly(subst_idx);
445 
446  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyBegin(bcvcorr[sbi]);
447  for (sbi=0; sbi<n_substances(); sbi++) VecAssemblyEnd(bcvcorr[sbi]);
448 
449  // we are calling set_boundary_conditions() after next_time() and
450  // we are using data from t() before, so we need to set corresponding bc time
452 }
453 
454 
455 //=============================================================================
456 // COMPUTE SOURCES
457 //=============================================================================
459 
460  //temporary variables
461  double csection, source, diag;
462 
463  //TODO: would it be possible to check the change in data for chosen substance? (may be in multifields?)
464 
465  //checking if the data were changed
466  if( (data_.sources_density.changed() )
467  || (data_.sources_conc.changed() )
468  || (data_.sources_sigma.changed() )
469  || (data_.cross_section.changed()))
470  {
471  START_TIMER("sources_reinit");
472  balance_->start_source_assembly(subst_idx);
473 
474  for ( DHCellAccessor dh_cell : dh_->own_range() )
475  {
476  ElementAccessor<3> elm = dh_cell.elm();
477  // we have currently zero order P_Disc FE
478  ASSERT_DBG(dh_cell.get_loc_dof_indices().size() == 1);
479  IntIdx local_p0_dof = dh_cell.get_loc_dof_indices()[0];
480 
481  arma::vec3 center = elm.centre();
482  csection = data_.cross_section.value(center, elm);
483 
484  // read for all substances
485  double max_cfl=0;
486  for (unsigned int sbi = 0; sbi < n_substances(); sbi++)
487  {
488  double src_sigma = data_.sources_sigma[sbi].value(center, elm);
489 
490  source = csection * (data_.sources_density[sbi].value(center, elm)
491  + src_sigma * data_.sources_conc[sbi].value(center, elm));
492  // addition to RHS
493  sources_corr[sbi][local_p0_dof] = source;
494  // addition to diagonal of the transport matrix
495  diag = src_sigma * csection;
496  tm_diag[sbi][local_p0_dof] = - diag;
497 
498  // compute maximal cfl condition over all substances
499  max_cfl = std::max(max_cfl, fabs(diag));
500 
501  balance_->add_source_values(sbi, elm.region().bulk_idx(), {local_p0_dof},
502  {- src_sigma * elm.measure() * csection},
503  {source * elm.measure()});
504  }
505 
506  cfl_source_[local_p0_dof] = max_cfl;
507  }
508 
509  balance_->finish_source_assembly(subst_idx);
510 
511  END_TIMER("sources_reinit");
512  }
513 }
514 
515 
516 
518 {
520 
523  std::stringstream ss; // print warning message with table of uninitialized fields
524  if ( FieldCommon::print_message_table(ss, "convection transport") ) {
525  WarningOut() << ss.str();
526  }
527 
530 
531  START_TIMER("Convection balance zero time step");
532 
536 
537  // write initial condition
538  output_data();
539 }
540 
541 
543 {
544  ASSERT_PTR(dh_).error( "Null DOF handler object.\n" );
545  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
546 
547  START_TIMER("data reinit");
548 
549  bool cfl_changed = false;
550 
551  // if FLOW or DATA changed ---------------------> recompute transport matrix
552  if (changed_)
553  {
556  cfl_changed = true;
557  DebugOut() << "CFL changed - flow.\n";
558  }
559 
561  {
563  cfl_changed = true;
564  DebugOut() << "CFL changed - mass matrix.\n";
565  }
566 
567  // if DATA changed ---------------------> recompute concentration sources (rhs and matrix diagonal)
570  {
572  is_src_term_scaled = false;
573  cfl_changed = true;
574  DebugOut() << "CFL changed - source.\n";
575  }
576 
577  // now resolve the CFL condition
578  if(cfl_changed)
579  {
580  // find maximum of sum of contribution from flow and sources: MAX(vcfl_flow_ + vcfl_source_)
581  Vec cfl;
582  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(),PETSC_DETERMINE, &cfl);
583  VecWAXPY(cfl, 1.0, vcfl_flow_, vcfl_source_);
584  VecMaxPointwiseDivide(cfl,mass_diag, &cfl_max_step);
585  // get a reciprocal value as a time constraint
587  DebugOut().fmt("CFL constraint (transport): {}\n", cfl_max_step);
588  }
589 
590  // although it does not influence CFL, compute BC so the full system is assembled
591  if ( changed_
592  || data_.porosity.changed()
594  || data_.bc_conc.changed() )
595  {
597  is_bc_term_scaled = false;
598  }
599 
600  END_TIMER("data reinit");
601 
602  // return time constraint
603  time_constraint = cfl_max_step;
604  changed_ = false;
605  return cfl_changed;
606 }
607 
609 
610  START_TIMER("convection-one step");
611 
612  // proceed to next time - which we are about to compute
613  // explicit scheme looks one step back and uses data from previous time
614  // (data time set previously in assess_time_constraint())
615  time_->next_time();
616 
617  double dt_new = time_->dt(), // current time step we are about to compute
618  dt_scaled = dt_new / time_->last_dt(); // scaling ratio to previous time step
619 
620  START_TIMER("time step rescaling");
621 
622  // if FLOW or DATA or BC or DT changed ---------------------> rescale boundary condition
624  {
625  DebugOut() << "BC - rescale dt.\n";
626  //choose between fresh scaling with new dt or rescaling to a new dt
627  double dt = (!is_bc_term_scaled) ? dt_new : dt_scaled;
628  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
629  VecScale(bcvcorr[sbi], dt);
630  is_bc_term_scaled = true;
631  }
632 
633 
634  // if DATA or TIME STEP changed -----------------------> rescale source term
636  DebugOut() << "SRC - rescale dt.\n";
637  //choose between fresh scaling with new dt or rescaling to a new dt
638  double dt = (!is_src_term_scaled) ? dt_new : dt_scaled;
639  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
640  {
641  VecScale(v_sources_corr[sbi], dt);
642  VecScale(v_tm_diag[sbi], dt);
643  }
644  is_src_term_scaled = true;
645  }
646 
647  // if DATA or TIME STEP changed -----------------------> rescale transport matrix
649  DebugOut() << "TM - rescale dt.\n";
650  //choose between fresh scaling with new dt or rescaling to a new dt
651  double dt = (!is_convection_matrix_scaled) ? dt_new : dt_scaled;
652 
653  MatScale(tm, dt);
655  }
656 
657  END_TIMER("time step rescaling");
658 
659 
660  data_.set_time(time_->step(), LimitSide::right); // set to the last computed time
662  {
663  VecCopy(mass_diag, vpmass_diag);
665  } else is_mass_diag_changed = false;
666 
667 
668  // Compute new concentrations for every substance.
669 
670  for (unsigned int sbi = 0; sbi < n_substances(); sbi++) {
671  // one step in MOBILE phase
672  START_TIMER("mat mult");
673 
674  // tm_diag is a diagonal part of transport matrix, which depends on substance data (sources_sigma)
675  // Wwe need keep transport matrix independent of substance, therefore we keep this diagonal part
676  // separately in a vector tm_diag.
677  // Computation: first, we compute this diagonal addition D*pconc and save it temporaly into RHS
678 
679  // RHS = D*pconc, where D is diagonal matrix represented by a vector
680  VecPointwiseMult(vcumulative_corr[sbi], v_tm_diag[sbi], vconc[sbi]); //w = x.*y
681 
682  // Then we add boundary terms ans other source terms into RHS.
683  // RHS = 1.0 * bcvcorr + 1.0 * v_sources_corr + 1.0 * rhs
684  VecAXPBYPCZ(vcumulative_corr[sbi], 1.0, 1.0, 1.0, bcvcorr[sbi], v_sources_corr[sbi]); //z = ax + by + cz
685 
686  // Then we set the new previous concentration.
687  VecCopy(vconc[sbi], vpconc[sbi]); // pconc = conc
688  // And finally proceed with transport matrix multiplication.
689  if (is_mass_diag_changed) {
690  VecPointwiseMult(vconc[sbi], vconc[sbi], vpmass_diag); // vconc*=vpmass_diag
691  MatMultAdd(tm, vpconc[sbi], vconc[sbi], vconc[sbi]); // vconc+=tm*vpconc
692  VecAXPY(vconc[sbi], 1, vcumulative_corr[sbi]); // vconc+=vcumulative_corr
693  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
694  } else {
695  MatMultAdd(tm, vpconc[sbi], vcumulative_corr[sbi], vconc[sbi]); // vconc =tm*vpconc+vcumulative_corr
696  VecPointwiseDivide(vconc[sbi], vconc[sbi], mass_diag); // vconc/=mass_diag
697  VecAXPY(vconc[sbi], 1, vpconc[sbi]); // vconc+=vpconc
698  }
699 
700  END_TIMER("mat mult");
701  }
702 
703  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
704  balance_->calculate_cumulative(sbi, vpconc[sbi]);
705 
706  END_TIMER("convection-one step");
707 }
708 
709 
710 void ConvectionTransport::set_target_time(double target_time)
711 {
712 
713  //sets target_mark_type (it is fixed) to be met in next_time()
714  time_->marks().add(TimeMark(target_time, target_mark_type));
715 
716  // This is done every time TOS calls update_solution.
717  // If CFL condition is changed, time fixation will change later from TOS.
718 
719  // Set the same constraint as was set last time.
720 
721  // TODO: fix this hack, remove this method completely, leaving just the first line at the calling point
722  // in TransportOperatorSplitting::update_solution()
723  // doing this directly leads to choose of large time step violationg CFL condition
724  if (cfl_max_step > time_->dt()*1e-10)
725  time_->set_upper_constraint(cfl_max_step, "CFL condition used from previous step.");
726 
727  // fixing convection time governor till next target_mark_type (got from TOS or other)
728  // may have marks for data changes
730 }
731 
732 
734 {
735  ElementAccessor<3> elm;
736 
737  VecZeroEntries(mass_diag);
738 
739  balance_->start_mass_assembly(subst_idx);
740 
741  for ( DHCellAccessor dh_cell : dh_->own_range() ) {
742  ElementAccessor<3> elm = dh_cell.elm();
743  // we have currently zero order P_Disc FE
744  ASSERT_DBG(dh_cell.get_loc_dof_indices().size() == 1);
745  IntIdx local_p0_dof = dh_cell.get_loc_dof_indices()[0];
746 
747  double csection = data_.cross_section.value(elm.centre(), elm);
748  //double por_m = data_.porosity.value(elm.centre(), elm->element_accessor());
749  double por_m = data_.water_content.value(elm.centre(), elm);
750 
751  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
752  balance_->add_mass_values(subst_idx[sbi], dh_cell, {local_p0_dof}, {csection*por_m*elm.measure()}, 0);
753 
754  VecSetValue(mass_diag, dh_->get_local_to_global_map()[local_p0_dof], csection*por_m, INSERT_VALUES);
755  }
756 
757  balance_->finish_mass_assembly(subst_idx);
758 
759  VecAssemblyBegin(mass_diag);
760  VecAssemblyEnd(mass_diag);
761 
762  is_mass_diag_changed = true;
763 }
764 
765 
766 //=============================================================================
767 // CREATE TRANSPORT MATRIX
768 //=============================================================================
770 
771  START_TIMER("convection_matrix_assembly");
772 
773  ElementAccessor<3> el2;
774  ElementAccessor<3> elm;
775  int j;
776  LongIdx new_j, new_i;
777  double aij, aii;
778 
779  MatZeroEntries(tm);
780 
781  double flux, flux2, edg_flux;
782 
783  aii = 0.0;
784 
785  unsigned int loc_el = 0;
786  for ( DHCellAccessor dh_cell : dh_->own_range() ) {
787  new_i = row_4_el[ dh_cell.elm_idx() ];
788  elm = dh_cell.elm();
789  for( DHCellSide cell_side : dh_cell.side_range() ) {
790  flux = this->side_flux(cell_side);
791  if (! cell_side.side().is_boundary()) {
792  edg_flux = 0;
793  for( DHCellSide edge_side : cell_side.edge_sides() ) {
794  el2 = edge_side.element();
795  flux2 = this->side_flux(edge_side);
796  if ( flux2 > 0) edg_flux+= flux2;
797  }
798  for( DHCellSide edge_side : cell_side.edge_sides() )
799  if (edge_side != cell_side) {
800  j = edge_side.element().idx();
801  new_j = row_4_el[j];
802 
803  el2 = edge_side.element();
804  flux2 = this->side_flux(edge_side);
805  if ( flux2 > 0.0 && flux <0.0)
806  aij = -(flux * flux2 / ( edg_flux * dh_cell.elm().measure() ) );
807  else aij =0;
808  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
809  }
810  }
811  if (flux > 0.0)
812  aii -= (flux / dh_cell.elm().measure() );
813  } // end same dim //ELEMENT_SIDES
814 
815  for( DHCellSide neighb_side : dh_cell.neighb_sides() ) // dh_cell lower dim
816  {
817  ASSERT( neighb_side.elem_idx() != dh_cell.elm_idx() ).error("Elm. same\n");
818  new_j = row_4_el[ neighb_side.elem_idx() ];
819  el2 = neighb_side.element();
820  flux = this->side_flux(neighb_side);
821 
822  // volume source - out-flow from higher dimension
823  if (flux > 0.0) aij = flux / dh_cell.elm().measure();
824  else aij=0;
825  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
826  // out flow from higher dim. already accounted
827 
828  // volume drain - in-flow to higher dimension
829  if (flux < 0.0) {
830  aii -= (-flux) / dh_cell.elm().measure(); // diagonal drain
831  aij = (-flux) / neighb_side.element().measure();
832  } else aij=0;
833  MatSetValue(tm, new_j, new_i, aij, INSERT_VALUES);
834  }
835 
836  MatSetValue(tm, new_i, new_i, aii, INSERT_VALUES);
837 
838  cfl_flow_[loc_el++] = fabs(aii);
839  aii = 0.0;
840  }
841 
842  MatAssemblyBegin(tm, MAT_FINAL_ASSEMBLY);
843  MatAssemblyEnd(tm, MAT_FINAL_ASSEMBLY);
844 
846  END_TIMER("convection_matrix_assembly");
847 
849 }
850 
851 
853  return conc;
854 }
855 
856 void ConvectionTransport::get_par_info(LongIdx * &el_4_loc_out, Distribution * &el_distribution_out){
857  el_4_loc_out = this->el_4_loc;
858  el_distribution_out = this->el_ds;
859  return;
860 }
861 
862 //int *ConvectionTransport::get_el_4_loc(){
863 // return el_4_loc;
864 //}
865 
867  return row_4_el;
868 }
869 
870 
871 
873 
875  data_.output_fields.output(time().step());
876 
877  START_TIMER("TOS-balance");
878  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
879  balance_->calculate_instant(sbi, vconc[sbi]);
880  balance_->output();
881  END_TIMER("TOS-balance");
882 }
883 
884 void ConvectionTransport::set_balance_object(std::shared_ptr<Balance> balance)
885 {
886  balance_ = balance;
887  subst_idx = balance_->add_quantities(substances_.names());
888 }
889 
891 {
892  if (cell_side.dim()==3) return calculate_side_flux<3>(cell_side);
893  else if (cell_side.dim()==2) return calculate_side_flux<2>(cell_side);
894  else return calculate_side_flux<1>(cell_side);
895 }
896 
897 template<unsigned int dim>
899 {
900  ASSERT_EQ(cell_side.dim(), dim).error("Element dimension mismatch!");
901 
902  feo_.fe_values(dim).reinit(cell_side.side());
903  auto vel = velocity_field_ptr_->value(cell_side.centre(), cell_side.element());
904  double side_flux = arma::dot(vel, feo_.fe_values(dim).normal_vector(0)) * feo_.fe_values(dim).JxW(0);
905  return side_flux;
906 }
TimeGovernor & time()
Definition: equation.hh:151
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:214
unsigned int size() const
get global size
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:313
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:608
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:338
int tlevel() const
double calculate_side_flux(const DHCellSide &cell)
Calculate flux on side of given element specified by dimension.
Definition: transport.cc:898
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:315
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:241
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:733
void next_time()
Proceed to the next time according to current estimated time step.
void initialize() override
Definition: transport.cc:182
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:362
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:320
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:301
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:339
double * cfl_flow_
Definition: transport.h:325
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:395
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:316
LongIdx * el_4_loc
Definition: transport.h:366
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
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > output_field_ptr
Fields correspond with out_conc.
Definition: transport.h:357
double side_flux(const DHCellSide &cell_side)
Wrapper of following method, call side_flux with correct template parameter.
Definition: transport.cc:890
const RegionDB & region_db() const
Definition: mesh.h:143
void set_balance_object(std::shared_ptr< Balance > balance) override
Definition: transport.cc:884
#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:344
const TimeStep & step(int index=-1) const
virtual ~ConvectionTransport()
Definition: transport.cc:254
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:517
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:360
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.
std::vector< VectorMPI > out_conc
Definition: transport.h:352
SubstanceList substances_
Transported substances.
Definition: transport.h:370
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:866
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:852
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:377
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:380
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:342
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:769
Support classes for parallel programing.
void alloc_transport_structs_mpi()
Definition: transport.cc:346
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:310
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:193
double * cfl_source_
Definition: transport.h:325
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)
bool changed_
Indicator of change in velocity field.
Definition: transport.h:386
#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:458
bool evaluate_time_constraint(double &time_constraint) override
Definition: transport.cc:542
FiniteElement< 1 > * fe1_
Definition: transport.h:102
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > velocity_field_ptr_
Pointer to velocity field given from Flow equation.
Definition: transport.h:383
double dt() const
Vec vcfl_flow_
Parallel vector for flow contribution to CFL condition.
Definition: transport.h:323
void set_target_time(double target_time) override
Definition: transport.cc:710
double ** tm_diag
Definition: transport.h:333
virtual void output_data() override
Write computed fields.
Definition: transport.cc:872
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:856
Distributed sparse graphs, partitioning.
#define ASSERT_DBG(expr)
Vec vcfl_source_
Parallel vector for source term contribution to CFL condition.
Definition: transport.h:323
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:350
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:389
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:367
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:227
LongIdx * row_4_el
Definition: transport.h:365
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:321
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.