Flow123d  jenkins-Flow123d-windows32-release-multijob-28
transport.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  *
17  * You should have received a copy of the GNU General Public License along with this program; if not,
18  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
19  *
20  *
21  * $Id$
22  * $Revision$
23  * $LastChangedBy$
24  * $LastChangedDate$
25  *
26  * @file
27  * @ingroup transport
28  * @brief Transport
29  *
30  *
31  */
32 
33 
34 #include <memory>
35 
36 #include "system/system.hh"
37 #include "system/sys_profiler.hh"
38 
39 #include "mesh/mesh.h"
40 #include "mesh/partitioning.hh"
41 #include "transport/transport.h"
42 
43 #include "io/output.h"
44 
45 #include "la/distribution.hh"
46 
47 #include "la/sparse_graph.hh"
48 #include <iostream>
49 #include <iomanip>
50 #include <string>
51 
52 // TODO: move partitioning into mesh_ and remove this include
53 #include "flow/darcy_flow_mh.hh"
54 #include "flow/old_bcd.hh"
55 #include "input/accessors.hh"
56 #include "input/input_type.hh"
57 
59 
61 #include "fields/field_values.hh"
63 #include "reaction/isotherm.hh" // SorptionType enum
64 
65 namespace IT = Input::Type;
66 
68  .add_value(none,"none","No sorption considered")
69  .add_value(Isotherm::linear,"linear","Linear isotherm described sorption considered.")
70  .add_value(Isotherm::freundlich,"freundlich","Freundlich isotherm described sorption considered")
71  .add_value(Isotherm::langmuir,"langmuir","Langmuir isotherm described sorption considered")
72  .close();
73 
74 
76  EqData().output_fields.make_output_field_selection("ConvectionTransport_Output")
77  .close();
78 
79 
81 {
82  ADD_FIELD(bc_conc, "Boundary conditions for concentrations.", "0.0");
83  ADD_FIELD(init_conc, "Initial concentrations.", "0.0");
84 
86 
87  output_fields += *this;
88  output_fields += conc_mobile.name("conc").units("M/L^3");
89 }
90 
91 
93 : TransportBase(init_mesh, in_rec)
94 {
95  this->eq_data_ = &data_;
96 
97  //mark type of the equation of convection transport (created in EquationBase constructor) and it is fixed
98  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"));
100 
102 
103  in_rec.val<Input::Array>("substances").copy_to(subst_names_);
104  n_subst_ = subst_names_.size();
105  INPUT_CHECK(n_subst_ >= 1 ,"Number of substances must be positive.\n");
106 
107  Input::Iterator<Input::Record> it = in_rec.find<Input::Record>("mass_balance");
108  if (it) mass_balance_ = new MassBalance(this, *it);
109 
111  data_.set_mesh(init_mesh);
112  data_.set_input_list( in_rec.val<Input::Array>("input_fields") );
114 
115 
116  sub_problem = 0;
117 
121  transport_matrix_time = -1.0; // or -infty
123  need_time_rescaling=true;
124 
125  // register output vectors
126  output_rec = in_rec.val<Input::Record>("output_stream");
130 
131  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
132  {
133  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
134  std::shared_ptr<FieldElementwise<3, FieldValue<3>::Scalar> > output_field_ptr(new FieldElementwise<3, FieldValue<3>::Scalar>(out_conc[MOBILE][sbi], n_subst_, mesh_->n_elements()));
135  data_.conc_mobile[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
136  }
141 
142 }
143 
144 
145 //=============================================================================
146 // MAKE TRANSPORT
147 //=============================================================================
149 
150  int * id_4_old = new int[mesh_->n_elements()];
151  int i = 0;
152  FOR_ELEMENTS(mesh_, ele) id_4_old[i++] = ele.index();
154  delete[] id_4_old;
155 
156  // TODO: make output of partitioning is usefull but makes outputs different
157  // on different number of processors, which breaks tests.
158  //
159  // Possible solution:
160  // - have flag in ini file to turn this output ON
161  // - possibility to have different ref_output for different num of proc.
162  // - or do not test such kind of output
163  //
164  //FOR_ELEMENTS(mesh_, ele) {
165  // ele->pid=el_ds->get_proc(row_4_el[ele.index()]);
166  //}
167 
168 }
169 
170 
171 
173 {
174  unsigned int sbi;
175 
176  if (mass_balance_ != NULL)
177  delete mass_balance_;
178 
179  //Destroy mpi vectors at first
180  VecDestroy(&v_sources_corr);
181  MatDestroy(&tm);
182 
183  VecDestroy(vconc);
184  VecDestroy(bcvcorr);
185  VecDestroy(vpconc);
186  VecDestroy(vcumulative_corr);
187  VecDestroy(vconc_out);
188 
189  for (sbi = 0; sbi < n_subst_; sbi++) {
190  //no mpi vectors
191  xfree(sources_density[sbi]);
192  xfree(sources_conc[sbi]);
193  xfree(sources_sigma[sbi]);
194  xfree(cumulative_corr[sbi]);
195  }
196 
197 
199 
204 
205  delete output_stream_;
206 
207 }
208 
209 
210 
211 
212 
214 {
215  FOR_ELEMENTS(mesh_, elem)
216  {
217  if (!el_ds->is_local(row_4_el[elem.index()])) continue;
218 
219  unsigned int index = row_4_el[elem.index()] - el_ds->begin();
220  ElementAccessor<3> ele_acc = mesh_->element_accessor(elem.index());
221  arma::vec value = data_.init_conc.value(elem->centre(), ele_acc);
222 
223  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
224  conc[MOBILE][sbi][index] = value(sbi);
225  }
226 
227 }
228 
229 //=============================================================================
230 // ALLOCATE OF TRANSPORT VARIABLES (ELEMENT & NODES)
231 //=============================================================================
233 
234  unsigned int i;
235  int sbi, n_subst, ph;
236  n_subst = n_subst_;
237 
238 
239  sources_corr = new double[el_ds->lsize()];
240  sources_density = (double**) xmalloc(n_subst * sizeof(double*));
241  sources_conc = (double**) xmalloc(n_subst * sizeof(double*));
242  sources_sigma = (double**) xmalloc(n_subst * sizeof(double*));
243 
244  cumulative_corr = (double**) xmalloc(n_subst * sizeof(double*));
245  for (sbi = 0; sbi < n_subst; sbi++) {
246  sources_density[sbi] = (double*) xmalloc(el_ds->lsize() * sizeof(double));
247  sources_conc[sbi] = (double*) xmalloc(el_ds->lsize() * sizeof(double));
248  sources_sigma[sbi] = (double*) xmalloc(el_ds->lsize() * sizeof(double));
249  cumulative_corr[sbi] = (double*) xmalloc(el_ds->lsize() * sizeof(double));
250  }
251 
252  conc = (double***) xmalloc(MAX_PHASES * sizeof(double**));
253  out_conc = (double***) xmalloc(MAX_PHASES * sizeof(double**));
254  for (ph = 0; ph < MAX_PHASES; ph++) {
255  if ((sub_problem & ph) == ph) {
256  conc[ph] = (double**) xmalloc(n_subst * sizeof(double*));
257  out_conc[ph] = (double**) xmalloc(n_subst * sizeof(double*));
258  for (sbi = 0; sbi < n_subst; sbi++) {
259  conc[ph][sbi] = (double*) xmalloc(el_ds->lsize() * sizeof(double));
260  out_conc[ph][sbi] = (double*) xmalloc(el_ds->size() * sizeof(double));
261  for (i = 0; i < el_ds->lsize(); i++) {
262  conc[ph][sbi][i] = 0.0;
263  }
264  for (i = 0; i < el_ds->size(); i++) {
265  out_conc[ph][sbi][i] = 0.0;
266  }
267  }
268  } else {
269  conc[ph] = NULL;
270  out_conc[ph] = NULL;
271  }
272  }
273 }
274 
275 //=============================================================================
276 // ALLOCATION OF TRANSPORT VECTORS (MPI)
277 //=============================================================================
279 
280  int sbi, n_subst, rank, np;
281  n_subst = n_subst_;
282 
283  MPI_Barrier(PETSC_COMM_WORLD);
284  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
285  MPI_Comm_size(PETSC_COMM_WORLD, &np);
286 
287  bcvcorr = (Vec*) xmalloc(n_subst * (sizeof(Vec)));
288  vconc = (Vec*) xmalloc(n_subst * (sizeof(Vec)));
289  vpconc = (Vec*) xmalloc(n_subst * (sizeof(Vec)));
290  vcumulative_corr = (Vec*) xmalloc(n_subst * (sizeof(Vec)));
291 
292 
293  vconc_out = (Vec*) xmalloc(n_subst * (sizeof(Vec))); // extend to all
294 
295 
296  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), PETSC_DECIDE,
298 
299  for (sbi = 0; sbi < n_subst; sbi++) {
300  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &bcvcorr[sbi]);
301  VecZeroEntries(bcvcorr[sbi]);
302  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(), conc[MOBILE][sbi],
303  &vconc[sbi]);
304 
305  VecCreateMPI(PETSC_COMM_WORLD, el_ds->lsize(), mesh_->n_elements(), &vpconc[sbi]);
306  VecZeroEntries(vconc[sbi]);
307  VecZeroEntries(vpconc[sbi]);
308 
309  // SOURCES
310  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, el_ds->lsize(), mesh_->n_elements(),
311  cumulative_corr[sbi],&vcumulative_corr[sbi]);
312 
313  VecCreateSeqWithArray(PETSC_COMM_SELF,1, mesh_->n_elements(), out_conc[MOBILE][sbi], &vconc_out[sbi]);
314 
315  VecZeroEntries(vcumulative_corr[sbi]);
316  VecZeroEntries(vconc_out[sbi]);
317  }
318 
319 
320  MatCreateAIJ(PETSC_COMM_WORLD, el_ds->lsize(), el_ds->lsize(), mesh_->n_elements(),
321  mesh_->n_elements(), 16, PETSC_NULL, 4, PETSC_NULL, &tm);
322 
323 }
324 
325 
327 {
329 
330  unsigned int sbi, loc_el;
331 
332  // Assembly bcvcorr vector
333  for(sbi=0; sbi < n_subst_; sbi++) VecZeroEntries(bcvcorr[sbi]);
334 
335 
336  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
337  elm = mesh_->element(el_4_loc[loc_el]);
338  if (elm->boundary_idx_ != NULL) {
339  unsigned int new_i = row_4_el[elm.index()];
340  double csection = data_.cross_section.value(elm->centre(), elm->element_accessor());
341  double por_m = data_.porosity.value(elm->centre(), elm->element_accessor());
342 
343  FOR_ELEMENT_SIDES(elm,si) {
344  Boundary *b = elm->side(si)->cond();
345  if (b != NULL) {
346  double flux = mh_dh->side_flux( *(elm->side(si)) );
347  if (flux < 0.0) {
348  double aij = -(flux / (elm->measure() * csection * por_m) );
349 
350  arma::vec value = data_.bc_conc.value( b->element()->centre(), b->element_accessor() );
351  for (sbi=0; sbi<n_subst_; sbi++)
352  VecSetValue(bcvcorr[sbi], new_i, value[sbi] * aij, ADD_VALUES);
353  }
354  }
355  }
356 
357  }
358  }
359 
360  for (sbi=0; sbi<n_subst_; sbi++) VecAssemblyBegin(bcvcorr[sbi]);
361  for (sbi=0; sbi<n_subst_; sbi++) VecAssemblyEnd(bcvcorr[sbi]);
362 
363 }
364 
365 
366 //=============================================================================
367 // COMPUTE SOURCES
368 //=============================================================================
370 
371  //temporary variables
372  unsigned int loc_el;
373  double conc_diff, csection;
374  ElementAccessor<3> ele_acc;
375  arma::vec3 p;
376 
377  //TODO: would it be possible to check the change in data for chosen substance? (may be in multifields?)
378 
379  //checking if the data were changed
380  if( (data_.sources_density.changed() )
381  || (data_.sources_conc.changed() )
382  || (data_.sources_sigma.changed() )
383  || (data_.cross_section.changed()))
384  {
385  START_TIMER("sources_reinit");
386  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++)
387  {
388  ele_acc = mesh_->element_accessor(el_4_loc[loc_el]);
389  p = ele_acc.centre();
390 
391  csection = data_.cross_section.value(p, ele_acc);
392 
393  //if(data_.sources_density.changed_during_set_time)
394  sources_density[sbi][loc_el] = data_.sources_density.value(p, ele_acc)(sbi)*csection;
395 
396  //if(data_.sources_conc.changed_during_set_time)
397  sources_conc[sbi][loc_el] = data_.sources_conc.value(p, ele_acc)(sbi);
398 
399  //if(data_.sources_sigma.changed_during_set_time)
400  sources_sigma[sbi][loc_el] = data_.sources_sigma.value(p, ele_acc)(sbi)*csection;
401  }
402  }
403 
404  //now computing source concentrations: density - sigma (source_conc - actual_conc)
405  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++)
406  {
407  conc_diff = sources_conc[sbi][loc_el] - conc[MOBILE][sbi][loc_el];
408  if ( conc_diff > 0.0)
409  sources_corr[loc_el] = ( sources_density[sbi][loc_el]
410  + conc_diff * sources_sigma[sbi][loc_el] )
411  * time_->dt();
412  else
413  sources_corr[loc_el] = sources_density[sbi][loc_el] * time_->dt();
414  }
415 }
416 
418 
419  //temporary variables
420  unsigned int loc_el;
421  double conc_diff, csection;
422  ElementAccessor<3> ele_acc;
423  arma::vec3 p;
424 
425  double *pconc;
426  VecGetArray(vpconc[sbi], &pconc);
427 
428  //TODO: would it be possible to check the change in data for chosen substance? (may be in multifields?)
429 
430  //checking if the data were changed
431  if( (data_.sources_density.changed() )
432  || (data_.sources_conc.changed() )
433  || (data_.sources_sigma.changed() )
434  || (data_.cross_section.changed()))
435  {
436  START_TIMER("sources_reinit");
437  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++)
438  {
439  ele_acc = mesh_->element_accessor(el_4_loc[loc_el]);
440  p = ele_acc.centre();
441 
442  csection = data_.cross_section.value(p, ele_acc);
443 
444  //if(data_.sources_density.changed_during_set_time)
445  sources_density[sbi][loc_el] = data_.sources_density.value(p, ele_acc)(sbi)*csection;
446 
447  //if(data_.sources_conc.changed_during_set_time)
448  sources_conc[sbi][loc_el] = data_.sources_conc.value(p, ele_acc)(sbi);
449 
450  //if(data_.sources_sigma.changed_during_set_time)
451  sources_sigma[sbi][loc_el] = data_.sources_sigma.value(p, ele_acc)(sbi)*csection;
452  }
453  }
454 
455  //now computing source concentrations: density - sigma (source_conc - actual_conc)
456  for (loc_el = 0; loc_el < el_ds->lsize(); loc_el++)
457  {
458  conc_diff = sources_conc[sbi][loc_el] - pconc[loc_el];
459  if ( conc_diff > 0.0)
460  sources_corr[loc_el] = sources_density[sbi][loc_el] + conc_diff * sources_sigma[sbi][loc_el];
461  else
462  sources_corr[loc_el] = sources_density[sbi][loc_el];
463  }
464 }
465 
466 
468 {
469  ASSERT_EQUAL(time_->tlevel(), 0);
470 
472  data_.set_time(*time_);
473 
475 
476 
477  // write initial condition
478  output_data();
479 }
480 
481 
482 
484 
485  START_TIMER("convection-one step");
486 
487  unsigned int sbi;
488 
489  START_TIMER("data reinit");
490  data_.set_time(*time_); // set to the last computed time
491 
492  ASSERT(mh_dh, "Null MH object.\n" );
493  // update matrix and sources if neccessary
494 
495 
499 
500  // need new fixation of the time step
501 
504 
506  // scale boundary sources
507  for (sbi=0; sbi<n_subst_; sbi++) VecScale(bcvcorr[sbi], time_->estimate_dt());
508 
509  need_time_rescaling = true;
510  } else {
511  // possibly read boundary conditions
512  if (data_.bc_conc.changed() ) {
514  // scale boundary sources
515  for (sbi=0; sbi<n_subst_; sbi++) VecScale(bcvcorr[sbi], time_->dt());
516  }
517  }
518 
519  if (need_time_rescaling) {
521  // rescale matrix
522  MatShift(tm, -1.0);
523  MatScale(tm, time_->estimate_dt()/time_->dt() );
524  MatShift(tm, 1.0);
525 
526  for (sbi=0; sbi<n_subst_; sbi++) VecScale(bcvcorr[sbi], time_->estimate_dt()/time_->dt());
527 
528  } else {
529  // scale fresh convection term matrix
530  MatScale(tm, time_->estimate_dt());
531  MatShift(tm, 1.0);
533 
534  }
535  need_time_rescaling = false;
536  }
537 
538  END_TIMER("data reinit");
539 
540 
541  // proceed to actually computed time
542  time_->next_time(); // explicit scheme use values from previous time and then set then new time
543 
544 
545  for (sbi = 0; sbi < n_subst_; sbi++) {
546  // one step in MOBILE phase
547 
548  START_TIMER("compute_concentration_sources");
549  //sources update
551 
552  //vcumulative_corr[sbi] = 1.0 * bcvcorr[sbi] + v_sources_corr;
553  VecWAXPY(vcumulative_corr[sbi],1.0,bcvcorr[sbi],v_sources_corr);
554  END_TIMER("compute_concentration_sources");
555 
556  START_TIMER("mat mult");
557  VecCopy(vconc[sbi], vpconc[sbi]); // pconc = conc
558  MatMultAdd(tm, vpconc[sbi], vcumulative_corr[sbi], vconc[sbi]); // conc=tm*pconc + bc
559  END_TIMER("mat mult");
560  }
561 
562  END_TIMER("convection-one step");
563 }
564 
565 
566 void ConvectionTransport::set_target_time(double target_time)
567 {
568 
569  //sets target_mark_type (it is fixed) to be met in next_time()
570  time_->marks().add(TimeMark(target_time, target_mark_type));
571 
572  // make new time step fixation, invalidate scaling
573  // same is done when matrix has changed in compute_one_step
575 
576  // fixing convection time governor till next target_mark_type (got from TOS or other)
577  // may have marks for data changes
579  need_time_rescaling = true;
580 
581 }
582 
583 
584 //=============================================================================
585 // CREATE TRANSPORT MATRIX
586 //=============================================================================
588 
589  START_TIMER("convection_matrix_assembly");
590 
593  struct Edge *edg;
594  unsigned int n;
595  int s, j, np, rank, new_j, new_i;
596  double max_sum, aij, aii;
597 
598  MatZeroEntries(tm);
599 
600  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
601  MPI_Comm_size(PETSC_COMM_WORLD, &np);
602 
603  double flux, flux2, edg_flux;
604 
605 
606  vector<double> edge_flow(mesh_->n_edges(),0.0);
607  for(unsigned int i=0; i < mesh_->n_edges() ; i++) { // calculate edge Qv
608  Edge &edg = mesh_->edges[i];
609  for( int s=0; s < edg.n_sides; s++) {
610  flux = mh_dh->side_flux( *(edg.side(s)) );
611  if ( flux > 0) edge_flow[i]+= flux;
612  }
613  }
614 
615  max_sum = 0.0;
616  aii = 0.0;
617 
618  for (unsigned int loc_el = 0; loc_el < el_ds->lsize(); loc_el++) {
619  elm = mesh_->element(el_4_loc[loc_el]);
620  new_i = row_4_el[elm.index()];
621 
622  double csection = data_.cross_section.value(elm->centre(), elm->element_accessor());
623  double por_m = data_.porosity.value(elm->centre(), elm->element_accessor());
624 
625  FOR_ELEMENT_SIDES(elm,si) {
626  // same dim
627  flux = mh_dh->side_flux( *(elm->side(si)) );
628  if (elm->side(si)->cond() == NULL) {
629  edg = elm->side(si)->edge();
630  edg_flux = edge_flow[ elm->side(si)->edge_idx() ];
631  FOR_EDGE_SIDES(edg,s)
632  // this test should also eliminate sides facing to lower dim. elements in comp. neighboring
633  // These edges on these sides should have just one side
634  if (edg->side(s) != elm->side(si)) {
635  j = ELEMENT_FULL_ITER(mesh_, edg->side(s)->element()).index();
636  new_j = row_4_el[j];
637 
638  flux2 = mh_dh->side_flux( *(edg->side(s)));
639  if ( flux2 > 0.0 && flux <0.0)
640  aij = -(flux * flux2 / ( edg_flux * elm->measure() * csection * por_m) );
641  else aij =0;
642  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
643  }
644  if (flux > 0.0)
645  aii -= (flux / (elm->measure() * csection * por_m) );
646  } else {
647  if (flux > 0.0)
648  aii -= (flux / (elm->measure() * csection * por_m) );
649  }
650  } // end same dim //ELEMENT_SIDES
651 
652  FOR_ELM_NEIGHS_VB(elm,n) // comp model
653  {
654  el2 = ELEMENT_FULL_ITER(mesh_, elm->neigh_vb[n]->side()->element() ); // higher dim. el.
655  ASSERT( el2 != elm, "Elm. same\n");
656  new_j = row_4_el[el2.index()];
657  flux = mh_dh->side_flux( *(elm->neigh_vb[n]->side()) );
658 
659  // volume source - out-flow from higher dimension
660  if (flux > 0.0) aij = flux / (elm->measure() * csection * por_m);
661  else aij=0;
662  MatSetValue(tm, new_i, new_j, aij, INSERT_VALUES);
663  // out flow from higher dim. already accounted
664 
665  // volume drain - in-flow to higher dimension
666  if (flux < 0.0) {
667  aii -= (-flux) / (elm->measure() * csection * por_m); // diagonal drain
668  aij = (-flux) / (el2->measure() *
669  data_.cross_section.value(el2->centre(), el2->element_accessor()) *
670  data_.porosity.value(el2->centre(), el2->element_accessor()));
671  } else aij=0;
672  MatSetValue(tm, new_j, new_i, aij, INSERT_VALUES);
673  }
674 
675  MatSetValue(tm, new_i, new_i, aii, INSERT_VALUES);
676 
677  if (fabs(aii) > max_sum)
678  max_sum = fabs(aii);
679  aii = 0.0;
680  } // END ELEMENTS
681 
682  double glob_max_sum;
683 
684  MPI_Allreduce(&max_sum,&glob_max_sum,1,MPI_DOUBLE,MPI_MAX,PETSC_COMM_WORLD);
685  cfl_max_step = 1 / glob_max_sum;
686 
687  MatAssemblyBegin(tm, MAT_FINAL_ASSEMBLY);
688  MatAssemblyEnd(tm, MAT_FINAL_ASSEMBLY);
689 
691  END_TIMER("convection_matrix_assembly");
692 
694 }
695 
696 
697 
698 
699 
700 //=============================================================================
701 // OUTPUT VECTOR GATHER
702 //=============================================================================
704 
705  unsigned int sbi;
706  IS is;
707 
708  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el, PETSC_COPY_VALUES, &is); //WithArray
709  VecScatterCreate(vconc[0], is, vconc_out[0], PETSC_NULL, &vconc_out_scatter);
710  for (sbi = 0; sbi < n_subst_; sbi++) {
711  VecScatterBegin(vconc_out_scatter, vconc[sbi], vconc_out[sbi], INSERT_VALUES, SCATTER_FORWARD);
712  VecScatterEnd(vconc_out_scatter, vconc[sbi], vconc_out[sbi], INSERT_VALUES, SCATTER_FORWARD);
713  }
714  VecScatterDestroy(&(vconc_out_scatter));
715  ISDestroy(&(is));
716 }
717 
718 
720  return conc;
721 }
722 
723 void ConvectionTransport::get_par_info(int * &el_4_loc_out, Distribution * &el_distribution_out){
724  el_4_loc_out = this->el_4_loc;
725  el_distribution_out = this->el_ds;
726  return;
727 }
728 
730  return el_4_loc;
731 }
732 
734  return row_4_el;
735 }
736 
737 
738 
739 void ConvectionTransport::calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance)
740 {
741  double mass_flux[n_substances()];
742  double *pconc[n_substances()];
743 
744  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
745  VecGetArray(vpconc[sbi], &pconc[sbi]);
746 
747  FOR_BOUNDARIES(mesh_, bcd) {
748 
749  // !! there can be more sides per one boundary
750  int index = row_4_el[bcd->side()->element().index()];
751  if (!el_ds->is_local(index)) continue;
752  int loc_index = index-el_ds->begin();
753 
754  double por_m = data_.porosity.value(bcd->side()->element()->centre(), bcd->side()->element()->element_accessor() );
755  double water_flux = mh_dh->side_flux(*(bcd->side()));
756  if (water_flux < 0) {
757  arma::vec bc_conc = data_.bc_conc.value( bcd->element()->centre(), bcd->element_accessor() );
758  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
759  mass_flux[sbi] = water_flux*bc_conc[sbi]*por_m;
760  } else {
761  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
762  mass_flux[sbi] = water_flux*pconc[sbi][loc_index]*por_m;
763  }
764 
765  Region r = bcd->region();
766  if (! r.is_valid()) xprintf(Msg, "Invalid region, ele % d, edg: % d\n", bcd->bc_ele_idx_, bcd->edge_idx_);
767  unsigned int bc_region_idx = r.boundary_idx();
768 
769  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
770  {
771  bcd_balance[sbi][bc_region_idx] += mass_flux[sbi];
772 
773  if (water_flux > 0) bcd_plus_balance[sbi][bc_region_idx] += mass_flux[sbi];
774  else bcd_minus_balance[sbi][bc_region_idx] += mass_flux[sbi];
775  }
776  }
777 
778 }
779 
781 {
782  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
783  {
785  double *sources = sources_corr;
786 
787  FOR_ELEMENTS(mesh_,elem)
788  {
789  int index = row_4_el[elem.index()];
790  if (!el_ds->is_local(index)) continue;
791  ElementAccessor<3> ele_acc = elem->element_accessor();
792  double por_m = data_.porosity.value(elem->centre(), ele_acc);
793  double csection = data_.cross_section.value(elem->centre(), ele_acc);
794  int loc_index = index - el_ds->begin();
795 
796  //temporary fix - mass balance works only when no reaction term is present
797  double sum_sol_phases = conc[0][sbi][loc_index];
798 
799  mass[sbi][ele_acc.region().bulk_idx()] += por_m*csection*sum_sol_phases*elem->measure();
800  src_balance[sbi][ele_acc.region().bulk_idx()] += sources[loc_index]*elem->measure();
801  }
802  }
803 }
804 
805 
806 
808 
809  if (time_->is_current( time_->marks().type_output() )) {
811 
814 
815 
817 
818  }
819 }
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:194
unsigned int size() const
get global size
void set_input_list(Input::Array input_list)
Definition: field_set.hh:164
void compute_concentration_sources(unsigned int sbi)
Definition: transport.cc:369
FieldSet * eq_data_
Definition: equation.hh:237
void set_limit_side(LimitSide side)
Definition: field_set.hh:171
double *** get_concentration_matrix()
Definition: transport.cc:719
double time_changed() const
double * sources_corr
Definition: transport.h:232
void init(const vector< string > &names)
Definition: field.impl.hh:487
Header: The functions for all outputs.
void update_solution() override
Definition: transport.cc:483
Accessor to input data conforming to declared Array.
Definition: accessors.hh:521
double end_time() const
End time.
Definition: system.hh:75
double transport_matrix_time
Definition: transport.h:253
int tlevel() const
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
void id_maps(int n_ids, int *id_4_old, Distribution *&new_ds, int *&id_4_loc, int *&new_4_id)
double *** conc
Concentrations for phase, substance, element.
Definition: transport.h:258
#define FOR_EDGE_SIDES(i, j)
Definition: edges.h:53
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:359
void alloc_transport_vectors()
Definition: transport.cc:232
double ** sources_sigma
Definition: transport.h:237
BCField< 3, FieldValue< 3 >::Vector > bc_conc
Definition: transport.h:101
FieldBasePtr(* read_field_descriptor_hook)(Input::Record rec, const FieldCommon &field)
Definition: field.hh:71
double fix_dt_until_mark()
Fixing time step until fixed time mark.
MultiField< 3, FieldValue< 3 >::Scalar > conc_mobile
Calculated concentrations in the mobile zone.
Definition: transport.h:106
void next_time()
Proceed to the next time according to current estimated time step.
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:241
TimeMark::Type type_output()
Definition: time_marks.hh:195
void set_initial_condition()
Definition: transport.cc:213
???
void set_boundary_conditions()
Definition: transport.cc:326
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:365
const MH_DofHandler * mh_dh
MassBalance * mass_balance_
object for calculation and writing the mass balance to file.
#define ELEMENT_FULL_ITER(_mesh_, i)
Definition: mesh.h:367
Definition: mesh.h:108
Field< 3, FieldValue< 3 >::Vector > sources_conc
Iterator< Ret > find(const string &key) const
#define FOR_ELEMENT_SIDES(i, j)
Definition: elements.h:151
int index() const
Definition: sys_vector.hh:88
bool is_current(const TimeMark::Type &mask) const
int n_sides
Definition: edges.h:48
Field< 3, FieldValue< 3 >::Scalar > porosity
Mobile porosity.
Definition: edges.h:38
const RegionDB & region_db() const
Definition: mesh.h:148
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:733
virtual ~ConvectionTransport()
Definition: transport.cc:172
double t() const
void get_par_info(int *&el_4_loc, Distribution *&el_ds)
Definition: transport.cc:723
#define ADD_FIELD(name,...)
Definition: field_set.hh:260
Basic time management functionality for unsteady (and steady) solvers (class Equation).
Partitioning * get_part()
Definition: mesh.cc:163
void zero_time_step() override
Definition: transport.cc:467
static TimeMarks & marks()
void output_vector_gather()
Definition: transport.cc:703
Basic time management class.
FieldCommon & units(const string &units)
Set basic units of the field.
Specification of transport model interface.
double estimate_dt() const
Estimate choice of next time step according to actual setting of constraints.
OutputTime * output_stream_
Definition: transport.h:272
static FieldBaseVector trans_conc_hook(Input::Record rec, const FieldCommon &field)
Definition: old_bcd.hh:117
ConvectionTransport(Mesh &init_mesh, const Input::Record &in_rec)
Definition: transport.cc:92
unsigned int boundary_idx() const
Returns index of the region in the boundary set.
Definition: region.hh:78
void compute_concentration_sources_for_mass_balance(unsigned int sbi)
Definition: transport.cc:417
void add(const TimeMark &mark)
Definition: time_marks.cc:85
#define MOBILE
unsigned int n_elements() const
Definition: mesh.h:137
#define ASSERT(...)
Definition: global_defs.h:120
bool is_local(unsigned int idx) const
identify local index
static OutputTime * create_output_stream(const Input::Record &in_rec)
This method write all registered data to output streams.
Definition: output.cc:171
double *** out_conc
Definition: transport.h:267
static Input::Type::Selection sorption_type_selection
Definition: transport.h:87
#define ASSERT_EQUAL(a, b)
Definition: global_defs.h:135
Input::Record output_rec
Record with output specification.
Definition: transport.h:270
void calc_fluxes(vector< vector< double > > &bcd_balance, vector< vector< double > > &bcd_plus_balance, vector< vector< double > > &bcd_minus_balance)
Definition: transport.cc:739
TimeMark::Type equation_fixed_mark_type() const
unsigned int begin(int proc) const
get starting local index
Selection & add_value(const int value, const std::string &key, const std::string &description="")
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:104
#define START_TIMER(tag)
Starts a timer with specified tag.
void output(double time)
Write computed fields to file.
Mesh * mesh_
Definition: equation.hh:227
Element * element()
Definition: boundaries.cc:47
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:399
void mark_input_times(TimeMark::Type mark_type)
Definition: field_set.hh:201
unsigned int n_subst_
Number of transported substances.
void mark_output_times(const TimeGovernor &tg)
Definition: output.cc:331
#define MPI_Comm_size
Definition: mpi.h:235
#define MPI_DOUBLE
Definition: mpi.h:156
Vec * vconc
Concentration vectors for mobile phase.
Definition: transport.h:256
void * xmalloc(size_t size)
Memory allocation with checking.
Definition: system.cc:218
#define INPUT_CHECK(i,...)
Debugging macros.
Definition: global_defs.h:61
void set_time(const TimeGovernor &time)
Definition: field_set.hh:186
Region region() const
Definition: accessors.hh:88
void calc_elem_sources(vector< vector< double > > &mass, vector< vector< double > > &src_balance)
Definition: transport.cc:780
#define MPI_Comm_rank
Definition: mpi.h:236
std::vector< string > subst_names_
Names of transported substances.
FieldSet output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
Definition: transport.h:109
void create_transport_matrix_mpi()
Definition: transport.cc:587
Support classes for parallel programing.
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Definition: mpi.h:612
Field< 3, FieldValue< 3 >::Vector > sources_sigma
Concentration sources - Robin type, in_flux = sources_sigma * (sources_conc - mobile_conc) ...
void alloc_transport_structs_mpi()
Definition: transport.cc:278
void set_n_components(unsigned int n_comp)
Definition: field_set.hh:151
double side_flux(const Side &side) const
temporary replacement for DofHandler accessor, flux through given side
bool is_convection_matrix_scaled
Definition: transport.h:224
#define ELEMENT_FULL_ITER_NULL(_mesh_)
Definition: mesh.h:373
static Input::Type::Selection output_selection
Definition: transport.h:89
#define MAX_PHASES
Definition: transport.h:64
void output(OutputTime *stream)
Definition: field_set.cc:137
const Selection & close() const
ElementFullIter element() const
Definition: side_impl.hh:43
int set_upper_constraint(double upper)
Sets upper constraint for the next time step estimating.
Field< 3, FieldValue< 3 >::Vector > sources_density
Concentration sources - density of substance source, only positive part is used.
VecScatter vconc_out_scatter
Definition: transport.h:247
double ** sources_density
temporary arrays to store constant values of fields over time interval
Definition: transport.h:237
double dt() const
void add_admissible_field_names(const Input::Array &in_array, const Input::Type::Selection &in_sel)
Registers names of output fields that can be written using this stream.
Definition: output.cc:231
arma::vec3 centre() const
Definition: elements.cc:138
virtual void output_data() override
Write computed fields.
Definition: transport.cc:807
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:209
Distributed sparse graphs, partitioning.
#define FOR_BOUNDARIES(_mesh_, i)
Definition: mesh.h:377
void set_target_time(double target_time)
Definition: transport.cc:566
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:83
FieldCommon & name(const string &name)
Definition: field_common.hh:71
#define END_TIMER(tag)
Ends a timer with specified tag.
arma::vec::fixed< spacedim > centre() const
Definition: accessors.hh:81
#define FOR_ELM_NEIGHS_VB(i, j)
Definition: elements.h:152
#define MPI_MAX
Definition: mpi.h:197
double ** cumulative_corr
Definition: transport.h:264
mixed-hybrid model of linear Darcy flow, possibly unsteady.
bool changed() const
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:157
Class used for marking specified times at which some events occur.
Definition: time_marks.hh:44
unsigned int n_substances() override
Returns number of trnasported substances.
bool is_valid() const
Returns false if the region has undefined/invalid value.
Definition: region.hh:70
unsigned int n_edges() const
Definition: mesh.h:145
#define xfree(p)
Definition: system.hh:112
#define MPI_Barrier(comm)
Definition: mpi.h:531
Distribution * el_ds
Definition: transport.h:277
SideIter side(const unsigned int i) const
Definition: edges.h:43
void set_mesh(const Mesh &mesh) override
Definition: field.impl.hh:534
Template for classes storing finite set of named values.
Field< 3, FieldValue< 3 >::Vector > init_conc
Initial concentrations.
Definition: transport.h:104
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:55
void make_transport_partitioning()
Definition: transport.cc:148
friend class MassBalance
Definition: mass_balance.hh:62
TimeGovernor * time_
Definition: equation.hh:228
double ** sources_conc
Definition: transport.h:237
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:198
unsigned int lsize(int proc) const
get local size