Flow123d  build_with_4.0.3-6210fd3
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"
32 
33 #include "la/distribution.hh"
34 
35 #include "la/sparse_graph.hh"
36 #include <iostream>
37 #include <iomanip>
38 #include <string>
39 
40 #include "io/output_time.hh"
41 #include "tools/time_governor.hh"
42 #include "tools/mixed.hh"
43 #include "coupling/balance.hh"
44 #include "input/accessors.hh"
45 #include "input/input_type.hh"
46 
48 #include "fields/field_values.hh"
49 #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 
70  Input::register_class< ConvectionTransport, Mesh &, const Input::Record >(_equation_name) +
72 
74 {
75  return IT::Record(_equation_name, "Finite volume method, explicit in time, for advection only solute transport.")
78  .declare_key("input_fields", IT::Array(
79  EqFields().make_field_descriptor_type(_equation_name)),
81  "")
82  .declare_key("output",
83  EqFields().output_fields.make_output_type(_equation_name, ""),
84  IT::Default("{ \"fields\": [ \"conc\" ] }"),
85  "Specification of output fields and output times.")
86  .close();
87 }
88 
89 
91 {
92  *this += bc_conc.name("bc_conc")
93  .description("Boundary condition for concentration of substances.")
94  .input_default("0.0")
95  .units( UnitSI().kg().m(-3) );
96 
97  *this += init_conc.name("init_conc")
98  .description("Initial values for concentration of substances.")
99  .input_default("0.0")
100  .units( UnitSI().kg().m(-3) );
101 
102  output_fields += *this;
103  output_fields += conc_mobile.name("conc")
104  .units( UnitSI().kg().m(-3) )
106  .description("Concentration solution.");
107  output_fields += region_id.name("region_id")
110  .description("Region ids.");
111  output_fields += subdomain.name("subdomain")
114  .description("Subdomain ids of the domain decomposition.");
115 
116  this->add_coords_field();
117  this->set_default_fieldset();
118 }
119 
120 
121 /********************************************************************************
122  * ConvectionTransport
123  */
125 : ConcentrationTransportBase(init_mesh, in_rec),
126  input_rec(in_rec)
127 {
128  START_TIMER("ConvectionTransport");
129  eq_data_ = make_shared<EqData>();
130  eq_fields_ = make_shared<EqFields>();
131  this->eq_fieldset_ = eq_fields_;
132 
133  eq_data_->transport_matrix_time = -1.0; // or -infty
134  eq_data_->transport_bc_time = -1.0;
135  eq_data_->is_convection_matrix_scaled = false;
136  is_src_term_scaled = false;
137  is_bc_term_scaled = false;
138 
139  //initialization of DOF handler
140  MixedPtr<FE_P_disc> fe(0);
141  eq_data_->dh_ = make_shared<DOFHandlerMultiDim>(init_mesh);
142  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>( &init_mesh, fe);
143  eq_data_->dh_->distribute_dofs(ds);
144  vcumulative_corr = nullptr;
145 
146 }
147 
149 {
151 
153 
154  eq_fields_->set_components(eq_data_->substances_.names());
155  eq_fields_->set_input_list( input_rec.val<Input::Array>("input_fields"), *time_ );
156  eq_fields_->set_mesh(*mesh_);
157 
160 
161  Input::Array user_fields_arr;
162  if (input_rec.opt_val("user_fields", user_fields_arr)) {
163  this->init_user_fields(user_fields_arr, eq_fields_->output_fields);
164  }
165 
166  // register output vectors
167  eq_fields_->output_fields.set_components(eq_data_->substances_.names());
168  eq_fields_->output_fields.set_mesh(*mesh_);
169  eq_fields_->output_fields.output_type(OutputTime::ELEM_DATA);
170  eq_fields_->conc_mobile.setup_components();
173  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
174  {
175  // create shared pointer to a FieldFE and push this Field to output_field on all regions
176  eq_fields_->conc_mobile_fe[sbi] = create_field_fe< 3, FieldValue<3>::Scalar >(eq_data_->dh_);
177  eq_fields_->conc_mobile[sbi].set(eq_fields_->conc_mobile_fe[sbi], 0);
178  }
179  //output_stream_->add_admissible_field_names(input_rec.val<Input::Array>("output_fields"));
180  //output_stream_->mark_output_times(*time_);
181 
182  eq_fields_->output_fields.initialize(output_stream_, mesh_, input_rec.val<Input::Record>("output"), time() );
183  //cout << "Transport." << endl;
184  //cout << time().marks();
185 
186  balance_->allocate(mesh_->get_el_ds()->lsize(), 1);
187  eq_data_->balance_ = this->balance();
188  eq_data_->set_time_governor(this->time_);
189  eq_data_->max_edg_sides = max(this->mesh_->max_edge_sides(1), max(this->mesh_->max_edge_sides(2), this->mesh_->max_edge_sides(3)));
190 
196 }
197 
198 
200 {
201  return eq_fields_->conc_mobile_fe[sbi]->vec().petsc_vec();
202 }
203 
204 
205 
207 {
208  unsigned int sbi;
209 
210  if (vcumulative_corr) {
211  //Destroy mpi vectors at first
212  chkerr(MatDestroy(&eq_data_->tm));
213  chkerr(VecDestroy(&eq_data_->mass_diag));
214  chkerr(VecDestroy(&vpmass_diag));
215 
216  for (sbi = 0; sbi < n_substances(); sbi++) {
217  // mpi vectors
218  chkerr(VecDestroy(&vpconc[sbi]));
219  chkerr(VecDestroy(&eq_data_->bcvcorr[sbi]));
220  chkerr(VecDestroy(&vcumulative_corr[sbi]));
221  }
222 
223  // arrays of mpi vectors
224  delete vpconc;
225  delete eq_data_->bcvcorr;
226  delete vcumulative_corr;
227 
228  // assembly objects
229  delete mass_assembly_;
230  delete init_cond_assembly_;
232  delete matrix_mpi_assembly_;
233  }
234 }
235 
236 
237 
238 
239 
240 //=============================================================================
241 // ALLOCATE OF TRANSPORT VARIABLES (ELEMENT & NODES)
242 //=============================================================================
244 
245  eq_fields_->conc_mobile_fe.resize(this->n_substances());
246 }
247 
248 //=============================================================================
249 // ALLOCATION OF TRANSPORT VECTORS (MPI)
250 //=============================================================================
252 
253  int sbi, n_subst, lsize;
254  n_subst = n_substances();
255  lsize = mesh_->get_el_ds()->lsize();
256 
257  MPI_Barrier(PETSC_COMM_WORLD);
258 
259  vpconc = new Vec[n_subst];
260  eq_data_->bcvcorr = new Vec[n_subst];
261  vcumulative_corr = new Vec[n_subst];
262  eq_data_->tm_diag.reserve(eq_data_->n_substances());
263  eq_data_->corr_vec.reserve(eq_data_->n_substances());
264 
265 
266  for (sbi = 0; sbi < n_subst; sbi++) {
267  VecCreateMPI(PETSC_COMM_WORLD, lsize, mesh_->n_elements(), &eq_data_->bcvcorr[sbi]);
268  VecZeroEntries(eq_data_->bcvcorr[sbi]);
269 
270  VecCreateMPI(PETSC_COMM_WORLD, lsize, mesh_->n_elements(), &vpconc[sbi]);
271  VecZeroEntries(vpconc[sbi]);
272 
273  // SOURCES
274  VecCreateMPI(PETSC_COMM_WORLD, lsize, mesh_->n_elements(), &vcumulative_corr[sbi]);
275 
276  eq_data_->corr_vec.emplace_back(lsize, PETSC_COMM_WORLD);
277 
278  eq_data_->tm_diag.emplace_back(lsize, PETSC_COMM_WORLD);
279 
280  VecZeroEntries(vcumulative_corr[sbi]);
281  }
282 
283 
284  MatCreateAIJ(PETSC_COMM_WORLD, lsize, lsize, mesh_->n_elements(),
285  mesh_->n_elements(), 16, PETSC_NULL, 4, PETSC_NULL, &eq_data_->tm);
286 
287  VecCreateMPI(PETSC_COMM_WORLD, lsize, mesh_->n_elements(), &eq_data_->mass_diag);
288  VecCreateMPI(PETSC_COMM_WORLD, lsize, mesh_->n_elements(), &vpmass_diag);
289 
290  eq_data_->alloc_transport_structs_mpi(lsize);
291 }
292 
293 
295 {
296  ASSERT_EQ(time_->tlevel(), 0);
297 
298  eq_fields_->mark_input_times(*time_);
299  eq_fields_->set_time(time_->step(), LimitSide::right);
300  std::stringstream ss; // print warning message with table of uninitialized fields
301  if ( FieldCommon::print_message_table(ss, "convection transport") ) {
302  WarningOut() << ss.str();
303  }
304 
305  //set_initial_condition();
307  //create_mass_matrix();
309 
310  START_TIMER("Convection balance zero time step");
311 
312  START_TIMER("convection_matrix_assembly");
313  //create_transport_matrix_mpi();
315  END_TIMER("convection_matrix_assembly");
316  START_TIMER("sources_reinit_set_bc");
317  //conc_sources_bdr_conditions();
319  END_TIMER("sources_reinit_set_bc");
320 
321  // write initial condition
322  output_data();
323 }
324 
325 
327 {
328  ASSERT_PTR(eq_data_->dh_).error( "Null DOF handler object.\n" );
329  // read changed status before setting time
330  bool changed_flux = eq_fields_->flow_flux.changed();
331  eq_fields_->set_time(time_->step(), LimitSide::right); // set to the last computed time
332 
333  START_TIMER("data reinit");
334 
335  bool cfl_changed = false;
336 
337  // if FLOW or DATA changed ---------------------> recompute transport matrix
338  if (changed_flux)
339  {
340  START_TIMER("convection_matrix_assembly");
341  //create_transport_matrix_mpi();
343  END_TIMER("convection_matrix_assembly");
344  eq_data_->is_convection_matrix_scaled=false;
345  cfl_changed = true;
346  DebugOut() << "CFL changed - flow.\n";
347  }
348 
349  if (eq_data_->is_mass_diag_changed)
350  {
351  //create_mass_matrix();
353  cfl_changed = true;
354  DebugOut() << "CFL changed - mass matrix.\n";
355  }
356 
357  // if DATA changed ---------------------> recompute concentration sources (rhs and matrix diagonal)
358  if( eq_fields_->sources_density.changed() || eq_fields_->sources_conc.changed() || eq_fields_->sources_sigma.changed()
359  || eq_fields_->cross_section.changed() || eq_fields_->flow_flux.changed() || eq_fields_->porosity.changed()
360  || eq_fields_->water_content.changed() || eq_fields_->bc_conc.changed() )
361  {
362  START_TIMER("sources_reinit_set_bc");
363  //conc_sources_bdr_conditions();
365  END_TIMER("sources_reinit_set_bc");
366  if( eq_data_->sources_changed_ ) {
367  is_src_term_scaled = false;
368  cfl_changed = true;
369  }
370  if (eq_fields_->flow_flux.changed() || eq_fields_->porosity.changed()
371  || eq_fields_->water_content.changed() || eq_fields_->bc_conc.changed() ) is_bc_term_scaled = false;
372  DebugOut() << "CFL changed - source.\n";
373  }
374 
375  // now resolve the CFL condition
376  if(cfl_changed)
377  {
378  // find maximum of sum of contribution from flow and sources: MAX(vcfl_flow_ + vcfl_source_)
379  Vec cfl;
380  VecCreateMPI(PETSC_COMM_WORLD, mesh_->get_el_ds()->lsize(), PETSC_DETERMINE, &cfl);
381  VecWAXPY(cfl, 1.0, eq_data_->cfl_flow_.petsc_vec(), eq_data_->cfl_source_.petsc_vec());
382  VecMaxPointwiseDivide(cfl, eq_data_->mass_diag, &cfl_max_step);
383  // get a reciprocal value as a time constraint
385  DebugOut().fmt("CFL constraint (transport): {}\n", cfl_max_step);
386  }
387 
388  END_TIMER("data reinit");
389 
390  // return time constraint
391  time_constraint = cfl_max_step;
392  return cfl_changed;
393 }
394 
396 
397  START_TIMER("convection-one step");
398 
399  // proceed to next time - which we are about to compute
400  // explicit scheme looks one step back and uses data from previous time
401  // (data time set previously in assess_time_constraint())
402  time_->next_time();
403 
404  double dt_new = time_->dt(), // current time step we are about to compute
405  dt_scaled = dt_new / time_->last_dt(); // scaling ratio to previous time step
406 
407  START_TIMER("time step rescaling");
408 
409  // if FLOW or DATA or BC or DT changed ---------------------> rescale boundary condition
411  {
412  DebugOut() << "BC - rescale dt.\n";
413  //choose between fresh scaling with new dt or rescaling to a new dt
414  double dt = (!is_bc_term_scaled) ? dt_new : dt_scaled;
415  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
416  VecScale(eq_data_->bcvcorr[sbi], dt);
417  is_bc_term_scaled = true;
418  }
419 
420 
421  // if DATA or TIME STEP changed -----------------------> rescale source term
423  DebugOut() << "SRC - rescale dt.\n";
424  //choose between fresh scaling with new dt or rescaling to a new dt
425  double dt = (!is_src_term_scaled) ? dt_new : dt_scaled;
426  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
427  {
428  VecScale(eq_data_->corr_vec[sbi].petsc_vec(), dt);
429  VecScale(eq_data_->tm_diag[sbi].petsc_vec(), dt);
430  }
431  is_src_term_scaled = true;
432  }
433 
434  // if DATA or TIME STEP changed -----------------------> rescale transport matrix
435  if ( !eq_data_->is_convection_matrix_scaled || time_->is_changed_dt()) {
436  DebugOut() << "TM - rescale dt.\n";
437  //choose between fresh scaling with new dt or rescaling to a new dt
438  double dt = (!eq_data_->is_convection_matrix_scaled) ? dt_new : dt_scaled;
439 
440  MatScale(eq_data_->tm, dt);
441  eq_data_->is_convection_matrix_scaled = true;
442  }
443 
444  END_TIMER("time step rescaling");
445 
446 
447  eq_fields_->set_time(time_->step(), LimitSide::right); // set to the last computed time
448  if (eq_fields_->cross_section.changed() || eq_fields_->water_content.changed() || eq_fields_->porosity.changed())
449  {
450  VecCopy(eq_data_->mass_diag, vpmass_diag);
451  //create_mass_matrix();
453  } else eq_data_->is_mass_diag_changed = false;
454 
455 
456  // Compute new concentrations for every substance.
457 
458  for (unsigned int sbi = 0; sbi < n_substances(); sbi++) {
459  // one step in MOBILE phase
460  START_TIMER("mat mult");
461  Vec vconc = eq_fields_->conc_mobile_fe[sbi]->vec().petsc_vec();
462 
463  // tm_diag is a diagonal part of transport matrix, which depends on substance data (sources_sigma)
464  // Wwe need keep transport matrix independent of substance, therefore we keep this diagonal part
465  // separately in a vector tm_diag.
466  // Computation: first, we compute this diagonal addition D*pconc and save it temporaly into RHS
467 
468  // RHS = D*pconc, where D is diagonal matrix represented by a vector
469  VecPointwiseMult(vcumulative_corr[sbi], eq_data_->tm_diag[sbi].petsc_vec(), vconc); //w = x.*y
470 
471  // Then we add boundary terms ans other source terms into RHS.
472  // RHS = 1.0 * bcvcorr + 1.0 * corr_vec + 1.0 * rhs
473  VecAXPBYPCZ(vcumulative_corr[sbi], 1.0, 1.0, 1.0, eq_data_->bcvcorr[sbi], eq_data_->corr_vec[sbi].petsc_vec()); //z = ax + by + cz
474 
475  // Then we set the new previous concentration.
476  VecCopy(vconc, vpconc[sbi]); // pconc = conc
477  // And finally proceed with transport matrix multiplication.
478  if (eq_data_->is_mass_diag_changed) {
479  VecPointwiseMult(vconc, vconc, vpmass_diag); // vconc*=vpmass_diag
480  MatMultAdd(eq_data_->tm, vpconc[sbi], vconc, vconc); // vconc+=tm*vpconc
481  VecAXPY(vconc, 1, vcumulative_corr[sbi]); // vconc+=vcumulative_corr
482  VecPointwiseDivide(vconc, vconc, eq_data_->mass_diag); // vconc/=mass_diag
483  } else {
484  MatMultAdd(eq_data_->tm, vpconc[sbi], vcumulative_corr[sbi], vconc); // vconc =tm*vpconc+vcumulative_corr
485  VecPointwiseDivide(vconc, vconc, eq_data_->mass_diag); // vconc/=mass_diag
486  VecAXPY(vconc, 1, vpconc[sbi]); // vconc+=vpconc
487  }
488 
489  END_TIMER("mat mult");
490  }
491 
492  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
493  balance_->calculate_cumulative(sbi, vpconc[sbi]);
494 
495  END_TIMER("convection-one step");
496 }
497 
498 
499 void ConvectionTransport::set_target_time(double target_time)
500 {
501 
502  //sets target_mark_type (it is fixed) to be met in next_time()
503  time_->marks().add(TimeMark(target_time, target_mark_type));
504 
505  // This is done every time TOS calls update_solution.
506  // If CFL condition is changed, time fixation will change later from TOS.
507 
508  // Set the same constraint as was set last time.
509 
510  // TODO: fix this hack, remove this method completely, leaving just the first line at the calling point
511  // in TransportOperatorSplitting::update_solution()
512  // doing this directly leads to choose of large time step violationg CFL condition
513  if (cfl_max_step > time_->dt()*1e-10)
514  time_->set_upper_constraint(cfl_max_step, "CFL condition used from previous step.");
515 
516  // fixing convection time governor till next target_mark_type (got from TOS or other)
517  // may have marks for data changes
519 }
520 
521 
523  return eq_fields_->conc_mobile_fe;
524 }
525 
526 
528 
529  eq_fields_->output_fields.set_time(time().step(), LimitSide::right);
530  eq_fields_->output_fields.output(time().step());
531 
532  START_TIMER("TOS-balance");
533  for (unsigned int sbi=0; sbi<n_substances(); ++sbi)
534  balance_->calculate_instant(sbi, eq_fields_->conc_mobile_fe[sbi]->vec().petsc_vec());
535  balance_->output();
536  END_TIMER("TOS-balance");
537 }
538 
539 void ConvectionTransport::set_balance_object(std::shared_ptr<Balance> balance)
540 {
541  balance_ = balance;
542  eq_data_->subst_idx = balance_->add_quantities(eq_data_->substances_.names());
543  eq_data_->balance_ = this->balance();
544 }
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Definition: asserts.hh:341
static Input::Type::Abstract & get_input_type()
Common specification of the input record for secondary equations.
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: transport.h:109
MultiField< 3, FieldValue< 3 >::Scalar > init_conc
Initial concentrations.
Definition: transport.h:106
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: transport.h:108
EquationOutput output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
Definition: transport.h:115
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_conc
Definition: transport.h:103
MultiField< 3, FieldValue< 3 >::Scalar > conc_mobile
Calculated concentrations in the mobile zone.
Definition: transport.h:111
FieldFEScalarVec & get_p0_interpolation() override
Getter for P0 interpolation by FieldFE.
Definition: transport.cc:522
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:271
void initialize() override
Definition: transport.cc:148
void set_target_time(double target_time) override
Definition: transport.cc:499
void update_solution() override
Definition: transport.cc:395
GenericAssembly< InitCondAssemblyConvection > * init_cond_assembly_
Definition: transport.h:323
void alloc_transport_structs_mpi()
Definition: transport.cc:251
std::shared_ptr< EqFields > eq_fields_
Definition: transport.h:291
void zero_time_step() override
Definition: transport.cc:294
GenericAssembly< MassAssemblyConvection > * mass_assembly_
general assembly objects, hold assembly objects of appropriate dimension
Definition: transport.h:322
void alloc_transport_vectors()
Definition: transport.cc:243
GenericAssembly< MatrixMpiAssemblyConvection > * matrix_mpi_assembly_
Definition: transport.h:325
ConvectionTransport(Mesh &init_mesh, const Input::Record in_rec)
Definition: transport.cc:124
const Input::Record input_rec
Record with input specification.
Definition: transport.h:316
virtual void output_data() override
Write computed fields.
Definition: transport.cc:527
Vec get_component_vec(unsigned int sbi) override
Return PETSc vector with solution for sbi-th substance.
Definition: transport.cc:199
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:318
void set_balance_object(std::shared_ptr< Balance > balance) override
Definition: transport.cc:539
static const Input::Type::Record & get_input_type()
Definition: transport.cc:73
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:304
GenericAssembly< ConcSourcesBdrAssemblyConvection > * conc_sources_bdr_assembly_
Definition: transport.h:324
bool evaluate_time_constraint(double &time_constraint) override
Definition: transport.cc:326
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:305
static const int registrar
Registrar of class to factory.
Definition: transport.h:286
std::shared_ptr< EqData > eq_data_
Definition: transport.h:292
virtual ~ConvectionTransport()
Definition: transport.cc:206
unsigned int lsize(int proc) const
get local size
std::shared_ptr< Balance > balance() const
Definition: equation.hh:189
std::shared_ptr< FieldSet > eq_fieldset_
Definition: equation.hh:249
void init_user_fields(Input::Array user_fields, FieldSet &output_fields)
Definition: equation.cc:84
TimeGovernor * time_
Definition: equation.hh:241
static Input::Type::Record & user_fields_template(std::string equation_name)
Template Record with common key user_fields for derived equations.
Definition: equation.cc:46
std::shared_ptr< Balance > balance_
object for calculation and writing the mass balance to file.
Definition: equation.hh:252
Mesh * mesh_
Definition: equation.hh:240
TimeGovernor & time()
Definition: equation.hh:151
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:95
FieldCommon & description(const string &description)
FieldCommon & flags(FieldFlag::Flags::Mask mask)
FieldCommon & name(const string &name)
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
FieldCommon & input_default(const string &input_default)
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:58
void add_coords_field()
Definition: field_set.cc:253
void set_default_fieldset()
Definition: field_set.hh:408
void assemble(std::shared_ptr< DOFHandlerMultiDim > dh) override
General assemble methods.
void set_min_edge_sides(unsigned int val)
static auto region_id(Mesh &mesh) -> IndexField
static auto subdomain(Mesh &mesh) -> IndexField
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
bool opt_val(const string &key, Ret &value) const
const Ret val(const string &key) const
Class for declaration of inputs sequences.
Definition: type_base.hh:339
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
Record type proxy class.
Definition: type_record.hh:182
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:602
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
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
Distribution * get_el_ds() const
Definition: mesh.h:119
unsigned int n_elements() const
Definition: mesh.h:111
unsigned int max_edge_sides(unsigned int dim) const
Definition: mesh.h:145
Definition: mesh.h:362
TimeMark::Type equation_fixed_mark_type() const
double dt() const
int set_upper_constraint(double upper, std::string message)
Sets upper constraint for the next time step estimating.
bool is_changed_dt() const
double end_time() const
End time.
double last_dt() const
double fix_dt_until_mark()
Fixing time step until fixed time mark.
const TimeStep & step(int index=-1) const
static TimeMarks & marks()
int tlevel() const
void next_time()
Proceed to the next time according to current estimated time step.
Class used for marking specified times at which some events occur.
Definition: time_marks.hh:45
TimeMark add(const TimeMark &mark)
Definition: time_marks.cc:81
Class for representation SI units of Fields.
Definition: unit_si.hh:40
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
Support classes for parallel programing.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:104
#define WarningOut()
Macro defining 'warning' record of log.
Definition: logger.hh:278
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
#define MPI_Barrier(comm)
Definition: mpi.h:531
Definitions of particular quadrature rules on simplices.
Implementation of range helper class.
Distributed sparse graphs, partitioning.
#define END_TIMER(tag)
Ends a timer with specified tag.
#define START_TIMER(tag)
Starts a timer with specified tag.
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:142
Basic time management class.
const string _equation_name
Definition: transport.cc:67