3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
transport.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
19#ifndef DUMUX_FVTRANSPORT_HH
20#define DUMUX_FVTRANSPORT_HH
21
22#include <dune/grid/common/gridenums.hh>
26#include <unordered_map>
27
33namespace Dumux
34{
36
50template<class TypeTag>
52{
54
56
57 enum
58 {
59 dim = GridView::dimension, dimWorld = GridView::dimensionworld
60 };
61
64
67
69
70 using Element = typename GridView::Traits::template Codim<0>::Entity;
71 using Intersection = typename GridView::Intersection;
72
73 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
74
75 struct LocalTimesteppingData
76 {
77 Dune::FieldVector<Scalar, 2*dim> faceFluxes;
78 Dune::FieldVector<Scalar, 2*dim> faceTargetDt;
79 Scalar dt;
80 LocalTimesteppingData():faceFluxes(0.0), faceTargetDt(0.0), dt(0)
81 {}
82 };
83
84protected:
86 EvalCflFluxFunction& evalCflFluxFunction()
87 {
88 return *evalCflFluxFunction_;
89 }
90
91 const EvalCflFluxFunction& evalCflFluxFunction() const
92 {
93 return *evalCflFluxFunction_;
94 }
96
97 void innerUpdate(TransportSolutionType& updateVec);
98
99public:
100
101 // Calculate the update vector.
102 void update(const Scalar t, Scalar& dt, TransportSolutionType& updateVec, bool impet = false);
103
104 void updateTransport(const Scalar t, Scalar& dt, TransportSolutionType& updateVec)
105 {
106 asImp_().updateMaterialLaws();
107 asImp_().update(t, dt, updateVec);
108 }
109
118 void getFlux(Scalar& update, const Intersection& intersection, CellData& cellDataI);
119
128 void getFluxOnBoundary(Scalar& update, const Intersection& intersection, CellData& cellDataI);
129
138 void getSource(Scalar& update, const Element& element, CellData& cellDataI);
139
142 {
143 evalCflFluxFunction_->initialize();
144 }
145
148
155 void getTransportedQuantity(TransportSolutionType& transportedQuantity);
156
157 template<class DataEntry>
158 bool inPhysicalRange(DataEntry& entry);
159
164 void setTransportedQuantity(TransportSolutionType& transportedQuantity);
165
170 void updateTransportedQuantity(TransportSolutionType& updateVec);
171
178 template<class MultiWriter>
179 void addOutputVtkFields(MultiWriter &writer)
180 {}
181
189 void serializeEntity(std::ostream &outstream, const Element &element)
190 {}
191
199 void deserializeEntity(std::istream &instream, const Element &element)
200 {}
201
203 {
204 return localTimeStepping_;
205 }
206
207
209
213 FVTransport(Problem& problem) :
214 problem_(problem), switchNormals_(getParam<bool>("Impet.SwitchNormals")),
215 subCFLFactor_(1.0), accumulatedDt_(0), dtThreshold_(1e-6)
216 {
217 evalCflFluxFunction_ = std::make_shared<EvalCflFluxFunction>(problem);
218
219 Scalar cFLFactor = getParam<Scalar>("Impet.CFLFactor");
220 using std::min;
221 subCFLFactor_ = min(getParam<Scalar>("Impet.SubCFLFactor"), cFLFactor);
222 verbosity_ = getParam<int>("TimeManager.SubTimestepVerbosity");
223
224 localTimeStepping_ = subCFLFactor_/cFLFactor < 1.0 - dtThreshold_;
225
226 if (localTimeStepping_)
227 std::cout<<"max CFL-Number of "<<cFLFactor<<", max Sub-CFL-Number of "<<subCFLFactor_<<": Enable local time-stepping!\n";
228 }
229
230private:
232 Implementation &asImp_()
233 { return *static_cast<Implementation *>(this); }
234
236 const Implementation &asImp_() const
237 { return *static_cast<const Implementation *>(this); }
238
239 void updatedTargetDt_(Scalar &dt);
240
241 void resetTimeStepData_()
242 {
243 timeStepData_.clear();
244 accumulatedDt_ = 0;
245 }
246
247 Problem& problem_;
248 bool switchNormals_;
249
250 std::shared_ptr<EvalCflFluxFunction> evalCflFluxFunction_;
251 std::vector<LocalTimesteppingData> timeStepData_;
252 bool localTimeStepping_;
253 Scalar subCFLFactor_;
254 Scalar accumulatedDt_;
255 const Scalar dtThreshold_;
256 int verbosity_;
257};
258
259
269template<class TypeTag>
270void FVTransport<TypeTag>::update(const Scalar t, Scalar& dt, TransportSolutionType& updateVec, bool impet)
271{
272 if (!impet)
273 {
274 asImp_().updateMaterialLaws();
275 }
276
277 unsigned int size = problem_.gridView().size(0);
278 if (localTimeStepping_)
279 {
280 if (timeStepData_.size() != size)
281 timeStepData_.resize(size);
282 }
283 // initialize dt very large
284 dt = std::numeric_limits<Scalar>::max();
285
286 // resize update vector and set to zero
287 updateVec.resize(size);
288 updateVec = 0.0;
289
290 // compute update vector
291 for (const auto& element : elements(problem_.gridView()))
292 {
293#if HAVE_MPI
294 if (element.partitionType() != Dune::InteriorEntity)
295 {
296 continue;
297 }
298#endif
299
300 // cell index
301 int globalIdxI = problem_.variables().index(element);
302
303 CellData& cellDataI = problem_.variables().cellData(globalIdxI);
304
305 Scalar update = 0;
306 evalCflFluxFunction().reset();
307
308 if (localTimeStepping_)
309 {
310 LocalTimesteppingData& localData = timeStepData_[globalIdxI];
311 for (int i = 0; i < 2*dim; i++)
312 {
313 if (localData.faceTargetDt[i] < accumulatedDt_ + dtThreshold_)
314 {
315 localData.faceFluxes[i] = 0.0;
316 }
317 }
318 }
319
320 // run through all intersections with neighbors and boundary
321 for (const auto& intersection : intersections(problem_.gridView(), element))
322 {
323 GlobalPosition unitOuterNormal = intersection.centerUnitOuterNormal();
324 if (switchNormals_)
325 unitOuterNormal *= -1.0;
326
327 int indexInInside = intersection.indexInInside();
328
329 // handle interior face
330 if (intersection.neighbor())
331 {
332 if (localTimeStepping_)
333 {
334 LocalTimesteppingData& localData = timeStepData_[globalIdxI];
335
336 if (localData.faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_)
337 {
338 asImp_().getFlux(localData.faceFluxes[indexInInside], intersection, cellDataI);
339 }
340 else
341 {
342 asImp_().getFlux(update, intersection, cellDataI);//only for time-stepping
343 }
344 }
345 else
346 {
347 //add flux to update
348 asImp_().getFlux(update, intersection, cellDataI);
349 }
350 } //end intersection with neighbor element
351 // handle boundary face
352 else if (intersection.boundary())
353 {
354 if (localTimeStepping_)
355 {
356 LocalTimesteppingData& localData = timeStepData_[globalIdxI];
357 if (localData.faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_)
358 {
359 asImp_().getFluxOnBoundary(localData.faceFluxes[indexInInside], intersection, cellDataI);
360 }
361 else
362 {
363 asImp_().getFluxOnBoundary(update, intersection, cellDataI);//only for time-stepping
364 }
365 }
366 else
367 {
368 //add boundary flux to update
369 asImp_().getFluxOnBoundary(update, intersection, cellDataI);
370 }
371 } //end boundary
372 } // end all intersections
373
374 if (localTimeStepping_)
375 {
376 LocalTimesteppingData& localData = timeStepData_[globalIdxI];
377 for (int i=0; i < 2*dim; i++)
378 {
379 updateVec[globalIdxI] += localData.faceFluxes[i];
380 }
381 }
382 else
383 {
384 //add flux update to global update vector
385 updateVec[globalIdxI] += update;
386 }
387
388 // std::cout<<"updateVec at "<<element.geometry().center()<<" : "<<updateVec[globalIdxI]<<"\n";
389
390 //add source to global update vector
391 Scalar source = 0.;
392 asImp_().getSource(source,element, cellDataI);
393 updateVec[globalIdxI] += source;
394
395 //calculate time step
396 using std::min;
397 if (localTimeStepping_)
398 {
399 Scalar dtCfl = evalCflFluxFunction().getDt(element);
400
401 timeStepData_[globalIdxI].dt = dtCfl;
402 dt = min(dt, dtCfl);
403 }
404 else
405 {
406 //calculate time step
407 dt = min(dt, evalCflFluxFunction().getDt(element));
408 }
409
410 //store update
411 cellDataI.setUpdate(updateVec[globalIdxI]);
412 } // end grid traversal
413
414
415#if HAVE_MPI
416 // communicate updated values
418 using ElementMapper = typename SolutionTypes::ElementMapper;
420 DataHandle dataHandle(problem_.elementMapper(), updateVec);
421 problem_.gridView().template communicate<DataHandle>(dataHandle,
422 Dune::InteriorBorder_All_Interface,
423 Dune::ForwardCommunication);
424
425 if (localTimeStepping_)
426 {
428
429 TimeDataHandle timeDataHandle(problem_.elementMapper(), timeStepData_);
430 problem_.gridView().template communicate<TimeDataHandle>(timeDataHandle,
431 Dune::InteriorBorder_All_Interface,
432 Dune::ForwardCommunication);
433 }
434
435 dt = problem_.gridView().comm().min(dt);
436#endif
437}
438
439template<class TypeTag>
441{
442 dt = std::numeric_limits<Scalar>::max();
443
444 // update target time-step-sizes
445 for (const auto& element : elements(problem_.gridView()))
446 {
447#if HAVE_MPI
448 if (element.partitionType() != Dune::InteriorEntity)
449 {
450 continue;
451 }
452#endif
453
454 // cell index
455 int globalIdxI = problem_.variables().index(element);
456
457 LocalTimesteppingData& localDataI = timeStepData_[globalIdxI];
458
459
460 using FaceDt = std::unordered_map<int, Scalar>;
461 FaceDt faceDt;
462
463 // run through all intersections with neighbors and boundary
464 for (const auto& intersection : intersections(problem_.gridView(), element))
465 {
466 int indexInInside = intersection.indexInInside();
467 using std::min;
468 if (intersection.neighbor())
469 {
470 auto neighbor = intersection.outside();
471 int globalIdxJ = problem_.variables().index(neighbor);
472
473 int levelI = element.level();
474 int levelJ = neighbor.level();
475
476 if (globalIdxI < globalIdxJ && levelI <= levelJ)
477 {
478 LocalTimesteppingData& localDataJ = timeStepData_[globalIdxJ];
479
480 int indexInOutside = intersection.indexInOutside();
481
482 if (localDataI.faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_
483 || localDataJ.faceTargetDt[indexInOutside] < accumulatedDt_ + dtThreshold_)
484 {
485 Scalar timeStep = min(localDataI.dt, localDataJ.dt);
486
487 if (levelI < levelJ)
488 {
489 typename FaceDt::iterator it = faceDt.find(indexInInside);
490 if (it != faceDt.end())
491 {
492 it->second = min(it->second, timeStep);
493 }
494 else
495 {
496 faceDt.insert(std::make_pair(indexInInside, timeStep));
497 }
498 }
499 else
500 {
501 localDataI.faceTargetDt[indexInInside] += subCFLFactor_ * timeStep;
502 localDataJ.faceTargetDt[indexInOutside] += subCFLFactor_ * timeStep;
503 }
504
505 dt = min(dt, timeStep);
506 }
507 }
508 }
509 else if (intersection.boundary())
510 {
511 if (localDataI.faceTargetDt[indexInInside] < accumulatedDt_ + dtThreshold_)
512 {
513 localDataI.faceTargetDt[indexInInside] += subCFLFactor_ * localDataI.dt;
514 dt =min(dt, subCFLFactor_ * localDataI.dt);
515 }
516 }
517 }
518 if (faceDt.size() > 0)
519 {
520 typename FaceDt::iterator it = faceDt.begin();
521 for (;it != faceDt.end();++it)
522 {
523 localDataI.faceTargetDt[it->first] += subCFLFactor_ * it->second;
524 }
525
526 for (const auto& intersection : intersections(problem_.gridView(), element))
527 {
528 if (intersection.neighbor())
529 {
530 int indexInInside = intersection.indexInInside();
531
532 it = faceDt.find(indexInInside);
533 if (it != faceDt.end())
534 {
535 int globalIdxJ = problem_.variables().index(intersection.outside());
536
537 LocalTimesteppingData& localDataJ = timeStepData_[globalIdxJ];
538
539 int indexInOutside = intersection.indexInOutside();
540
541 localDataJ.faceTargetDt[indexInOutside] += subCFLFactor_ * it->second;
542 }
543 }
544 }
545 }
546 }
547
548#if HAVE_MPI
549 // communicate updated values
550 using SolutionTypes = GetProp<TypeTag, Properties::SolutionTypes>;
551 using ElementMapper = typename SolutionTypes::ElementMapper;
552 using TimeDataHandle = VectorExchange<ElementMapper, std::vector<LocalTimesteppingData> >;
553
554 TimeDataHandle timeDataHandle(problem_.elementMapper(), timeStepData_);
555 problem_.gridView().template communicate<TimeDataHandle>(timeDataHandle,
556 Dune::InteriorBorder_All_Interface,
557 Dune::ForwardCommunication);
558
559 dt = problem_.gridView().comm().min(dt);
560#endif
561}
562
563template<class TypeTag>
564void FVTransport<TypeTag>::innerUpdate(TransportSolutionType& updateVec)
565{
566 if (localTimeStepping_)
567 {
568 Scalar realDt = problem_.timeManager().timeStepSize();
569
570 Scalar subDt = realDt;
571
572 updatedTargetDt_(subDt);
573
574 Scalar accumulatedDtOld = accumulatedDt_;
575 accumulatedDt_ += subDt;
576
577 Scalar t = problem_.timeManager().time();
578
579 if (accumulatedDt_ < realDt)
580 {
581 using std::min;
582 while(true)
583 {
584 Scalar dtCorrection = min(0.0, realDt - accumulatedDt_);
585 subDt += dtCorrection;
586
587 if (verbosity_ > 0)
588 std::cout<<" Sub-time-step size: "<<subDt<<"\n";
589
590 TransportSolutionType transportedQuantity;
591 asImp_().getTransportedQuantity(transportedQuantity);
592
593 bool stopTimeStep = false;
594 int size = transportedQuantity.size();
595 for (int i = 0; i < size; i++)
596 {
597 Scalar newVal = transportedQuantity[i] += updateVec[i] * subDt;
598 if (!asImp_().inPhysicalRange(newVal))
599 {
600 stopTimeStep = true;
601
602 break;
603 }
604 }
605
606#if HAVE_MPI
607 int rank = 0;
608 if (stopTimeStep)
609 rank = problem_.gridView().comm().rank();
610
611 rank = problem_.gridView().comm().max(rank);
612 problem_.gridView().comm().broadcast(&stopTimeStep,1,rank);
613#endif
614
615
616 if (stopTimeStep && accumulatedDtOld > dtThreshold_)
617 {
618 problem_.timeManager().setTimeStepSize(accumulatedDtOld);
619 break;
620 }
621 else
622 {
623 asImp_().setTransportedQuantity(transportedQuantity);
624 }
625
626
627 if (accumulatedDt_ >= realDt)
628 {
629 break;
630 }
631
632 problem_.model().updateTransport(t, subDt, updateVec);
633
634 updatedTargetDt_(subDt);
635
636 accumulatedDtOld = accumulatedDt_;
637 accumulatedDt_ += subDt;
638 }
639 }
640 else
641 {
642 asImp_().updateTransportedQuantity(updateVec, realDt);
643 }
644
645 resetTimeStepData_();
646 }
647}
648}
649#endif
Contains a class to exchange entries of a vector.
Specifies the properties for immiscible 2p transport.
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition: parameters.hh:365
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:79
The finite volume discretization of a transport equation.
Definition: transport.hh:52
void initialize()
Sets the initial solution .
Definition: transport.hh:141
void updateMaterialLaws()
Updates constitutive relations and stores them in the variable class.
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the primary transport variable.
Definition: transport.hh:199
void getTransportedQuantity(TransportSolutionType &transportedQuantity)
Writes the current values of the primary transport variable into the transportedQuantity-vector (come...
void addOutputVtkFields(MultiWriter &writer)
Adds transport output to the output file.
Definition: transport.hh:179
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the primary transport variable.
Definition: transport.hh:189
void getSource(Scalar &update, const Element &element, CellData &cellDataI)
Function which calculates the source update.
bool enableLocalTimeStepping()
Definition: transport.hh:202
void getFlux(Scalar &update, const Intersection &intersection, CellData &cellDataI)
Function which calculates the flux update.
void setTransportedQuantity(TransportSolutionType &transportedQuantity)
Writes the current values of the primary transport variable into the variable container.
void updateTransportedQuantity(TransportSolutionType &updateVec)
Updates the primary transport variable.
FVTransport(Problem &problem)
Constructs a FVTransport object.
Definition: transport.hh:213
void getFluxOnBoundary(Scalar &update, const Intersection &intersection, CellData &cellDataI)
Function which calculates the boundary flux update.
void updateTransport(const Scalar t, Scalar &dt, TransportSolutionType &updateVec)
Definition: transport.hh:104
void update(const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impet=false)
Calculate the update vector.
Definition: transport.hh:270
void innerUpdate(TransportSolutionType &updateVec)
Definition: transport.hh:564
bool inPhysicalRange(DataEntry &entry)
Base file for properties related to sequential models.