27#ifndef EWOMS_VTK_DISCRETE_FRACTURE_MODULE_HH
28#define EWOMS_VTK_DISCRETE_FRACTURE_MODULE_HH
30#include <dune/common/fvector.hh>
32#include <opm/material/common/Valgrind.hpp>
71template <
class TypeTag>
86 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
89 enum { dim = GridView::dimension };
90 enum { dimWorld = GridView::dimensionworld };
91 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
107 Parameters::Register<Parameters::VtkWriteFractureSaturations>
108 (
"Include the phase saturations in the VTK output files");
109 Parameters::Register<Parameters::VtkWriteFractureMobilities>
110 (
"Include the phase mobilities in the VTK output files");
111 Parameters::Register<Parameters::VtkWriteFractureRelativePermeabilities>
112 (
"Include the phase relative permeabilities in the "
114 Parameters::Register<Parameters::VtkWriteFracturePorosity>
115 (
"Include the porosity in the VTK output files");
116 Parameters::Register<Parameters::VtkWriteFractureIntrinsicPermeabilities>
117 (
"Include the intrinsic permeability in the VTK output files");
118 Parameters::Register<Parameters::VtkWriteFractureFilterVelocities>
119 (
"Include in the filter velocities of the phases in the VTK output files");
120 Parameters::Register<Parameters::VtkWriteFractureVolumeFraction>
121 (
"Add the fraction of the total volume which is "
122 "occupied by fractures in the VTK output");
131 if (saturationOutput_())
133 if (mobilityOutput_())
135 if (relativePermeabilityOutput_())
138 if (porosityOutput_())
140 if (intrinsicPermeabilityOutput_())
142 if (volumeFractionOutput_())
145 if (velocityOutput_()) {
146 size_t nDof = this->
simulator_.model().numGridDof();
147 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
148 fractureVelocity_[phaseIdx].resize(nDof);
149 for (
unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx) {
150 fractureVelocity_[phaseIdx][dofIdx].resize(dimWorld);
151 fractureVelocity_[phaseIdx][dofIdx] = 0.0;
164 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
168 const auto& fractureMapper = elemCtx.simulator().vanguard().fractureMapper();
170 for (
unsigned i = 0; i < elemCtx.numPrimaryDof(0); ++i) {
171 unsigned I = elemCtx.globalSpaceIndex(i, 0);
172 if (!fractureMapper.isFractureVertex(I))
175 const auto& intQuants = elemCtx.intensiveQuantities(i, 0);
176 const auto& fs = intQuants.fractureFluidState();
178 if (porosityOutput_()) {
179 Opm::Valgrind::CheckDefined(intQuants.fracturePorosity());
180 fracturePorosity_[I] = intQuants.fracturePorosity();
182 if (intrinsicPermeabilityOutput_()) {
183 const auto& K = intQuants.fractureIntrinsicPermeability();
184 fractureIntrinsicPermeability_[I] = K[0][0];
187 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
188 if (saturationOutput_()) {
189 Opm::Valgrind::CheckDefined(fs.saturation(phaseIdx));
190 fractureSaturation_[phaseIdx][I] = fs.saturation(phaseIdx);
192 if (mobilityOutput_()) {
193 Opm::Valgrind::CheckDefined(intQuants.fractureMobility(phaseIdx));
194 fractureMobility_[phaseIdx][I] = intQuants.fractureMobility(phaseIdx);
196 if (relativePermeabilityOutput_()) {
197 Opm::Valgrind::CheckDefined(intQuants.fractureRelativePermeability(phaseIdx));
198 fractureRelativePermeability_[phaseIdx][I] =
199 intQuants.fractureRelativePermeability(phaseIdx);
201 if (volumeFractionOutput_()) {
202 Opm::Valgrind::CheckDefined(intQuants.fractureVolume());
203 fractureVolumeFraction_[I] += intQuants.fractureVolume();
208 if (velocityOutput_()) {
210 for (
unsigned scvfIdx = 0; scvfIdx < elemCtx.numInteriorFaces(0); ++ scvfIdx) {
211 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, 0);
213 unsigned i = extQuants.interiorIndex();
214 unsigned I = elemCtx.globalSpaceIndex(i, 0);
216 unsigned j = extQuants.exteriorIndex();
217 unsigned J = elemCtx.globalSpaceIndex(j, 0);
219 if (!fractureMapper.isFractureEdge(I, J))
222 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
224 std::max<Scalar>(1e-16, std::abs(extQuants.fractureVolumeFlux(phaseIdx)));
225 Opm::Valgrind::CheckDefined(extQuants.extrusionFactor());
226 assert(extQuants.extrusionFactor() > 0);
227 weight *= extQuants.extrusionFactor();
229 Dune::FieldVector<Scalar, dim> v(extQuants.fractureFilterVelocity(phaseIdx));
232 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
233 fractureVelocity_[phaseIdx][I][dimIdx] += v[dimIdx];
234 fractureVelocity_[phaseIdx][J][dimIdx] += v[dimIdx];
237 fractureVelocityWeight_[phaseIdx][I] += weight;
238 fractureVelocityWeight_[phaseIdx][J] += weight;
254 if (saturationOutput_())
256 if (mobilityOutput_())
258 if (relativePermeabilityOutput_())
259 this->
commitPhaseBuffer_(baseWriter,
"fractureRelativePerm_%s", fractureRelativePermeability_);
261 if (porosityOutput_())
263 if (intrinsicPermeabilityOutput_())
264 this->
commitScalarBuffer_(baseWriter,
"fractureIntrinsicPerm", fractureIntrinsicPermeability_);
265 if (volumeFractionOutput_()) {
267 for (
unsigned I = 0; I < fractureVolumeFraction_.size(); ++I)
268 fractureVolumeFraction_[I] /= this->
simulator_.model().dofTotalVolume(I);
272 if (velocityOutput_()) {
273 size_t nDof = this->
simulator_.model().numGridDof();
275 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
278 for (
unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx)
279 fractureVelocity_[phaseIdx][dofIdx] /=
280 std::max<Scalar>(1e-20, fractureVelocityWeight_[phaseIdx][dofIdx]);
283 snprintf(name, 512,
"fractureFilterVelocity_%s", FluidSystem::phaseName(phaseIdx).data());
285 DiscBaseOutputModule::attachVectorDofData_(baseWriter, fractureVelocity_[phaseIdx], name);
291 static bool saturationOutput_()
293 static bool val = Parameters::Get<Parameters::VtkWriteFractureSaturations>();
297 static bool mobilityOutput_()
299 static bool val = Parameters::Get<Parameters::VtkWriteFractureMobilities>();
303 static bool relativePermeabilityOutput_()
305 static bool val = Parameters::Get<Parameters::VtkWriteFractureRelativePermeabilities>();
309 static bool porosityOutput_()
311 static bool val = Parameters::Get<Parameters::VtkWriteFracturePorosity>();
315 static bool intrinsicPermeabilityOutput_()
317 static bool val = Parameters::Get<Parameters::VtkWriteFractureIntrinsicPermeabilities>();
321 static bool volumeFractionOutput_()
323 static bool val = Parameters::Get<Parameters::VtkWriteFractureVolumeFraction>();
327 static bool velocityOutput_()
329 static bool val = Parameters::Get<Parameters::VtkWriteFractureFilterVelocities>();
333 PhaseBuffer fractureSaturation_;
334 PhaseBuffer fractureMobility_;
335 PhaseBuffer fractureRelativePermeability_;
337 ScalarBuffer fracturePorosity_;
338 ScalarBuffer fractureVolumeFraction_;
339 ScalarBuffer fractureIntrinsicPermeability_;
341 PhaseVectorBuffer fractureVelocity_;
342 PhaseBuffer fractureVelocityWeight_;
344 PhaseVectorBuffer potentialGradient_;
345 PhaseBuffer potentialWeight_;
The base class for writer modules.
Definition: baseoutputmodule.hh:67
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:85
void commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase-specific buffer to the result file.
Definition: baseoutputmodule.hh:345
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase-specific quantity.
Definition: baseoutputmodule.hh:201
const Simulator & simulator_
Definition: baseoutputmodule.hh:434
std::array< VectorBuffer, numPhases > PhaseVectorBuffer
Definition: baseoutputmodule.hh:94
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:90
void commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer containing scalar quantities to the result file.
Definition: baseoutputmodule.hh:244
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:156
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for quantities which make sense for all models which deal with discrete fractures i...
Definition: vtkdiscretefracturemodule.hh:73
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiscretefracturemodule.hh:129
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkdiscretefracturemodule.hh:105
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiscretefracturemodule.hh:247
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkdiscretefracturemodule.hh:162
VtkDiscreteFractureModule(const Simulator &simulator)
Definition: vtkdiscretefracturemodule.hh:98
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilnewtonmethodparameters.hh:31
Definition: blackoilboundaryratevector.hh:37
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:235
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: vtkdiscretefracturemodule.hh:52
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:52
Definition: vtkdiscretefracturemodule.hh:51
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:51
Definition: vtkdiscretefracturemodule.hh:48
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:48
Definition: vtkdiscretefracturemodule.hh:50
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:50
Definition: vtkdiscretefracturemodule.hh:49
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:49
Definition: vtkdiscretefracturemodule.hh:47
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:47
Definition: vtkdiscretefracturemodule.hh:53
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:53