Forked from
ogs / ogs
8245 commits behind the upstream repository.
-
Dmitri Naumov authoredDmitri Naumov authored
Code owners
Assign users and groups as approvers for specific file changes. Learn more.
CreateTwoPhaseFlowWithPPProcess.cpp 4.00 KiB
/**
* \file
* \copyright
* Copyright (c) 2012-2021, OpenGeoSys Community (http://www.opengeosys.org)
* Distributed under a Modified BSD License.
* See accompanying file LICENSE.txt or
* http://www.opengeosys.org/project/license
*
*/
#include "CreateTwoPhaseFlowWithPPProcess.h"
#include <cassert>
#include "MaterialLib/MPL/CreateMaterialSpatialDistributionMap.h"
#include "MaterialLib/MPL/MaterialSpatialDistributionMap.h"
#include "ParameterLib/Utils.h"
#include "ProcessLib/Output/CreateSecondaryVariables.h"
#include "ProcessLib/Utils/ProcessUtils.h"
#include "TwoPhaseFlowWithPPProcess.h"
#include "TwoPhaseFlowWithPPProcessData.h"
namespace ProcessLib
{
namespace TwoPhaseFlowWithPP
{
std::unique_ptr<Process> createTwoPhaseFlowWithPPProcess(
std::string name,
MeshLib::Mesh& mesh,
std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
std::vector<ProcessVariable> const& variables,
std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
unsigned const integration_order,
BaseLib::ConfigTree const& config,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves,
std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
{
//! \ogs_file_param{prj__processes__process__type}
config.checkConfigParameter("type", "TWOPHASE_FLOW_PP");
DBUG("Create TwoPhaseFlowProcess with PP model.");
//! \ogs_file_param{prj__processes__process__TWOPHASE_FLOW_PP__process_variables}
auto const pv_config = config.getConfigSubtree("process_variables");
auto per_process_variables = findProcessVariables(
variables, pv_config,
{//! \ogs_file_param_special{prj__processes__process__TWOPHASE_FLOW_PP__process_variables__gas_pressure}
"gas_pressure",
//! \ogs_file_param_special{prj__processes__process__TWOPHASE_FLOW_PP__process_variables__capillary_pressure}
"capillary_pressure"});
std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
process_variables;
process_variables.push_back(std::move(per_process_variables));
SecondaryVariableCollection secondary_variables;
ProcessLib::createSecondaryVariables(config, secondary_variables);
// Specific body force
std::vector<double> const b =
//! \ogs_file_param{prj__processes__process__TWOPHASE_FLOW_PP__specific_body_force}
config.getConfigParameter<std::vector<double>>("specific_body_force");
assert(!b.empty() && b.size() < 4);
Eigen::VectorXd specific_body_force(b.size());
bool const has_gravity = MathLib::toVector(b).norm() > 0;
if (has_gravity)
{
std::copy_n(b.data(), b.size(), specific_body_force.data());
}
//! \ogs_file_param{prj__processes__process__TWOPHASE_FLOW_PP__mass_lumping}
auto const mass_lumping = config.getConfigParameter<bool>("mass_lumping");
auto const& temperature = ParameterLib::findParameter<double>(
config,
//! \ogs_file_param_special{prj__processes__process__TWOPHASE_FLOW_PP__temperature}
"temperature", parameters, 1, &mesh);
auto const material_ids = materialIDs(mesh);
if (material_ids)
{
INFO("The twophase flow is in heterogeneous porous media.");
}
else
{
INFO("The twophase flow is in homogeneous porous media.");
}
auto media_map =
MaterialPropertyLib::createMaterialSpatialDistributionMap(media, mesh);
TwoPhaseFlowWithPPProcessData process_data{specific_body_force, has_gravity,
mass_lumping, temperature,
std::move(media_map)};
return std::make_unique<TwoPhaseFlowWithPPProcess>(
std::move(name), mesh, std::move(jacobian_assembler), parameters,
integration_order, std::move(process_variables),
std::move(process_data), std::move(secondary_variables), curves);
}
} // namespace TwoPhaseFlowWithPP
} // namespace ProcessLib