OGS
ProcessLib::ComponentTransport Namespace Reference

Classes

class  ComponentTransportLocalAssemblerInterface
class  ComponentTransportProcess
struct  ComponentTransportProcessData
struct  Field
struct  IntegrationPointData
struct  InterpolationPoint
class  LocalAssemblerData
struct  LookupTable

Functions

void checkMPLProperties (MeshLib::Mesh const &mesh, MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map)
std::unique_ptr< ProcesscreateComponentTransportProcess (std::string const &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::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media, std::unique_ptr< ChemistryLib::ChemicalSolverInterface > &&chemical_solver_interface)
std::unique_ptr< LookupTablecreateLookupTable (std::optional< std::string > const tabular_file, std::filesystem::path const &project_directory, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const &process_variables)
static void intersection (std::vector< std::size_t > &vec1, std::vector< std::size_t > const &vec2)

Function Documentation

◆ checkMPLProperties()

void ProcessLib::ComponentTransport::checkMPLProperties ( MeshLib::Mesh const & mesh,
MaterialPropertyLib::MaterialSpatialDistributionMap const & media_map )

Definition at line 24 of file CreateComponentTransportProcess.cpp.

27{
28 std::array const required_properties_medium = {
33
34 std::array const required_properties_liquid_phase = {
37
38 std::array const required_properties_components = {
42
43 for (auto const& element_id : mesh.getElements() | MeshLib::views::ids)
44 {
45 auto const& medium = *media_map.getMedium(element_id);
46 checkRequiredProperties(medium, required_properties_medium);
47
48 // check if liquid phase definition and the corresponding properties
49 // exist
50 auto const& liquid_phase =
52 checkRequiredProperties(liquid_phase, required_properties_liquid_phase);
53
54 // check if components and the corresponding properties exist
55 auto const number_of_components = liquid_phase.numberOfComponents();
56 for (std::size_t component_id = 0; component_id < number_of_components;
57 ++component_id)
58 {
59 if (!liquid_phase.hasComponent(component_id))
60 {
62 "The component {:d} in the AqueousLiquid phase isn't "
63 "specified.",
64 component_id);
65 }
66 auto const& component = liquid_phase.component(component_id);
67 checkRequiredProperties(component, required_properties_components);
68 }
69 }
70}
#define OGS_FATAL(...)
Definition Error.h:19
void checkRequiredProperties(Component const &c, std::span< PropertyType const > const required_properties)
Definition Component.cpp:51
@ longitudinal_dispersivity
used to compute the hydrodynamic dispersion tensor.
@ transversal_dispersivity
used to compute the hydrodynamic dispersion tensor.
@ retardation_factor
specify retardation factor used in component transport process.
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:218

References MaterialPropertyLib::AqueousLiquid, MaterialPropertyLib::decay_rate, MaterialPropertyLib::density, MeshLib::Mesh::getElements(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), MeshLib::views::ids, MaterialPropertyLib::longitudinal_dispersivity, OGS_FATAL, MaterialPropertyLib::permeability, MaterialPropertyLib::pore_diffusion, MaterialPropertyLib::porosity, MaterialPropertyLib::retardation_factor, MaterialPropertyLib::transversal_dispersivity, and MaterialPropertyLib::viscosity.

Referenced by createComponentTransportProcess().

◆ createComponentTransportProcess()

std::unique_ptr< Process > ProcessLib::ComponentTransport::createComponentTransportProcess ( std::string const & 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::vector< std::unique_ptr< MeshLib::Mesh > > const & meshes,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media,
std::unique_ptr< ChemistryLib::ChemicalSolverInterface > && chemical_solver_interface )
Input File Parameter
prj__processes__process__type
Input File Parameter
prj__processes__process__ComponentTransport__coupling_scheme

Process Variables

Input File Parameter
prj__processes__process__ComponentTransport__process_variables

Primary process variables as they appear in the global component vector.

Input File Parameter
prj__processes__process__ComponentTransport__process_variables__pressure
Input File Parameter
prj__processes__process__ComponentTransport__process_variables__concentration

Temperature is optional.

Input File Parameter
prj__processes__process__ComponentTransport__process_variables__temperature

Process Parameters

Input File Parameter
prj__processes__process__ComponentTransport__specific_body_force
Input File Parameter
prj__processes__process__ComponentTransport__non_advective_form
Input File Parameter
prj__processes__process__ComponentTransport__chemically_induced_porosity_change
Input File Parameter
prj__processes__process__ComponentTransport__temperature_field
Input File Parameter
prj__processes__process__ComponentTransport__tabular_file
Input File Parameter
prj__processes__process__ComponentTransport__aperture_size
Input File Parameter
prj__processes__process__ComponentTransport__aperture_size__parameter
Input File Parameter
prj__processes__process__linear
Input File Parameter
prj__processes__process__ComponentTransport__linear_solver_compute_only_upon_timestep_change
Input File Parameter
prj__processes__process__calculatesurfaceflux

Definition at line 72 of file CreateComponentTransportProcess.cpp.

84{
86 config.checkConfigParameter("type", "ComponentTransport");
87
88 DBUG("Create ComponentTransportProcess.");
89
90 auto const coupling_scheme =
92 config.getConfigParameter<std::string>("coupling_scheme",
93 "monolithic_scheme");
94 const bool use_monolithic_scheme = (coupling_scheme != "staggered");
95
97
99 auto const pv_config = config.getConfigSubtree("process_variables");
100
101 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
102 process_variables;
103
105 std::vector collected_process_variables = findProcessVariables(
106 variables, pv_config,
107 {
108 "pressure",
110 "concentration"});
111
113 std::vector const temperature_variable = findProcessVariables(
114 variables, pv_config,
116 "temperature", true /*temperature is optional*/);
117 bool const isothermal = temperature_variable.empty();
118 if (!isothermal)
119 {
120 assert(temperature_variable.size() == 1);
121 collected_process_variables.insert(
122 ++collected_process_variables.begin(), temperature_variable[0]);
123 }
124
125 // Check number of components for each process variable
126 auto it = std::find_if(
127 collected_process_variables.cbegin(),
128 collected_process_variables.cend(),
129 [](std::reference_wrapper<ProcessLib::ProcessVariable> const& pv)
130 { return pv.get().getNumberOfGlobalComponents() != 1; });
131
132 if (it != collected_process_variables.end())
133 {
134 OGS_FATAL(
135 "Number of components for process variable '{:s}' should be 1 "
136 "rather "
137 "than {:d}.",
138 it->get().getName(),
139 it->get().getNumberOfGlobalComponents());
140 }
141
142 // Allocate the collected process variables into a two-dimensional vector,
143 // depending on what scheme is adopted
144 if (use_monolithic_scheme) // monolithic scheme.
145 {
146 if (!isothermal)
147 {
148 OGS_FATAL(
149 "Currently, non-isothermal component transport process can "
150 "only be simulated in staggered scheme.");
151 }
152
153 process_variables.push_back(std::move(collected_process_variables));
154 }
155 else // staggered scheme.
156 {
157 std::vector<std::reference_wrapper<ProcessLib::ProcessVariable>>
158 per_process_variable;
159
160 if (!chemical_solver_interface)
161 {
162 for (auto& pv : collected_process_variables)
163 {
164 per_process_variable.emplace_back(pv);
165 process_variables.push_back(std::move(per_process_variable));
166 }
167 }
168 else
169 {
170 auto sort_by_component =
171 [&per_process_variable,
172 collected_process_variables](auto const& c_name)
173 {
174 auto pv = std::find_if(collected_process_variables.begin(),
175 collected_process_variables.end(),
176 [&c_name](auto const& v) -> bool
177 { return v.get().getName() == c_name; });
178
179 if (pv == collected_process_variables.end())
180 {
181 OGS_FATAL(
182 "Component {:s} given in "
183 "<chemical_system>/<solution>/"
184 "<components> is not found in specified "
185 "coupled processes (see "
186 "<process>/<process_variables>/"
187 "<concentration>).",
188 c_name);
189 }
190
191 per_process_variable.emplace_back(*pv);
192 return std::move(per_process_variable);
193 };
194
195 auto const components =
196 chemical_solver_interface->getComponentList();
197 // pressure
198 per_process_variable.emplace_back(collected_process_variables[0]);
199 process_variables.push_back(std::move(per_process_variable));
200 // temperature
201 if (!isothermal)
202 {
203 per_process_variable.emplace_back(
204 collected_process_variables[1]);
205 process_variables.push_back(std::move(per_process_variable));
206 }
207 // concentration
208 assert(components.size() + (isothermal ? 1 : 2) ==
209 collected_process_variables.size());
210 std::transform(components.begin(), components.end(),
211 std::back_inserter(process_variables),
212 sort_by_component);
213 }
214 }
215
217 std::vector<double> const b =
219 config.getConfigParameter<std::vector<double>>("specific_body_force");
220 assert(!b.empty() && b.size() < 4);
221 // Specific body force parameter.
222 Eigen::VectorXd specific_body_force(b.size());
223 int const mesh_space_dimension =
225 if (static_cast<int>(b.size()) != mesh_space_dimension)
226 {
227 OGS_FATAL(
228 "specific body force (gravity vector) has {:d} components, mesh "
229 "dimension is {:d}",
230 b.size(), mesh_space_dimension);
231 }
232 bool const has_gravity = MathLib::toVector(b).norm() > 0;
233 if (has_gravity)
234 {
235 std::copy_n(b.data(), b.size(), specific_body_force.data());
236 }
237
238 bool const non_advective_form =
240 config.getConfigParameter<bool>("non_advective_form", false);
241
242 bool chemically_induced_porosity_change =
244 config.getConfigParameter<bool>("chemically_induced_porosity_change",
245 false);
246
247 auto const temperature_field = ParameterLib::findOptionalTagParameter<
248 double>(
250 config, "temperature_field", parameters, 1, &mesh);
251 if (!isothermal && temperature_field != nullptr)
252 {
253 OGS_FATAL("Temperature field is set for non-isothermal setup.")
254 }
255
256 auto media_map =
258
259 auto lookup_table = ComponentTransport::createLookupTable(
261 config.getConfigParameterOptional<std::string>("tabular_file"),
262 config.projectDirectory(),
263 process_variables);
264
265 DBUG("Check the media properties of ComponentTransport process ...");
266 checkMPLProperties(mesh, media_map);
267 DBUG("Media properties verified.");
268
269 auto stabilizer = NumLib::createNumericalStabilization(mesh, config);
270
271 auto const* aperture_size_parameter = &ParameterLib::findParameter<double>(
273 auto const aperture_config =
275 config.getConfigSubtreeOptional("aperture_size");
276 if (aperture_config)
277 {
278 aperture_size_parameter = &ParameterLib::findParameter<double>(
280 *aperture_config, "parameter", parameters, 1);
281 }
282
283 auto const is_linear =
285 config.getConfigParameter("linear", false);
286
287 auto const ls_compute_only_upon_timestep_change =
289 config.getConfigParameter(
290 "linear_solver_compute_only_upon_timestep_change", false);
291
292 auto const rotation_matrices = MeshLib::getElementRotationMatrices(
293 mesh_space_dimension, mesh.getDimension(), mesh.getElements());
294 std::vector<Eigen::VectorXd> projected_specific_body_force_vectors;
295 projected_specific_body_force_vectors.reserve(rotation_matrices.size());
296
297 std::transform(rotation_matrices.begin(), rotation_matrices.end(),
298 std::back_inserter(projected_specific_body_force_vectors),
299 [&specific_body_force](const auto& R)
300 { return R * R.transpose() * specific_body_force; });
301
303 std::move(media_map),
304 has_gravity,
305 non_advective_form,
306 temperature_field,
307 chemically_induced_porosity_change,
308 chemical_solver_interface.get(),
309 std::move(lookup_table),
310 std::move(stabilizer),
311 projected_specific_body_force_vectors,
312 mesh_space_dimension,
313 *aperture_size_parameter,
314 isothermal,
315 NumLib::ShapeMatrixCache{integration_order}};
316
317 SecondaryVariableCollection secondary_variables;
318
319 ProcessLib::createSecondaryVariables(config, secondary_variables);
320
321 std::unique_ptr<ProcessLib::SurfaceFluxData> surfaceflux;
322 auto surfaceflux_config =
324 config.getConfigSubtreeOptional("calculatesurfaceflux");
325 if (surfaceflux_config)
326 {
328 *surfaceflux_config, meshes);
329 }
330
331 return std::make_unique<ComponentTransportProcess>(
332 std::move(name), mesh, std::move(jacobian_assembler), parameters,
333 integration_order, std::move(process_variables),
334 std::move(process_data), std::move(secondary_variables),
335 use_monolithic_scheme, std::move(surfaceflux),
336 std::move(chemical_solver_interface), is_linear,
337 ls_compute_only_upon_timestep_change);
338}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:98
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:101
unsigned getDimension() const
Definition Mesh.h:80
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name
Definition Process.h:40
MaterialSpatialDistributionMap createMaterialSpatialDistributionMap(std::map< int, std::shared_ptr< Medium > > const &media, MeshLib::Mesh const &mesh)
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
std::vector< Eigen::MatrixXd > getElementRotationMatrices(int const space_dimension, int const mesh_dimension, std::vector< Element * > const &elements)
Element rotation matrix computation.
int getSpaceDimension(std::vector< Node * > const &nodes)
Computes dimension of the embedding space containing the set of given points.
NumericalStabilization createNumericalStabilization(MeshLib::Mesh const &mesh, BaseLib::ConfigTree const &config)
Parameter< ParameterDataType > * findOptionalTagParameter(BaseLib::ConfigTree const &process_config, std::string const &tag, std::vector< std::unique_ptr< ParameterBase > > const &parameters, int const num_components, MeshLib::Mesh const *const mesh=nullptr)
OGS_NO_DANGLING Parameter< ParameterDataType > & findParameter(std::string const &parameter_name, std::vector< std::unique_ptr< ParameterBase > > const &parameters, int const num_components, MeshLib::Mesh const *const mesh=nullptr)
std::unique_ptr< LookupTable > createLookupTable(std::optional< std::string > const tabular_file, std::filesystem::path const &project_directory, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const &process_variables)
void checkMPLProperties(MeshLib::Mesh const &mesh, MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map)
std::vector< std::reference_wrapper< ProcessVariable > > findProcessVariables(std::vector< ProcessVariable > const &variables, BaseLib::ConfigTree const &pv_config, std::initializer_list< std::string > tags)
void createSecondaryVariables(BaseLib::ConfigTree const &config, SecondaryVariableCollection &secondary_variables)
static std::unique_ptr< ProcessLib::SurfaceFluxData > createSurfaceFluxData(BaseLib::ConfigTree const &calculatesurfaceflux_config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes)

References BaseLib::ConfigTree::checkConfigParameter(), checkMPLProperties(), ProcessLib::Process::constant_one_parameter_name, createLookupTable(), MaterialPropertyLib::createMaterialSpatialDistributionMap(), NumLib::createNumericalStabilization(), ProcessLib::createSecondaryVariables(), ProcessLib::SurfaceFluxData::createSurfaceFluxData(), DBUG(), ParameterLib::findOptionalTagParameter(), ParameterLib::findParameter(), ProcessLib::findProcessVariables(), BaseLib::ConfigTree::getConfigParameter(), BaseLib::ConfigTree::getConfigParameterOptional(), BaseLib::ConfigTree::getConfigSubtree(), BaseLib::ConfigTree::getConfigSubtreeOptional(), MeshLib::Mesh::getDimension(), MeshLib::getElementRotationMatrices(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::getSpaceDimension(), OGS_FATAL, BaseLib::ConfigTree::projectDirectory(), and MathLib::toVector().

Referenced by ProjectData::parseProcesses().

◆ createLookupTable()

std::unique_ptr< LookupTable > ProcessLib::ComponentTransport::createLookupTable ( std::optional< std::string > const tabular_file,
std::filesystem::path const & project_directory,
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & process_variables )

Definition at line 48 of file CreateLookupTable.cpp.

53{
54 if (!tabular_file)
55 {
56 return nullptr;
57 }
58
59 auto const path_to_tabular_file =
60 BaseLib::joinPaths(project_directory.string(), *tabular_file);
61
62 if (!BaseLib::IsFileExisting(path_to_tabular_file))
63 {
65 "Not found the tabular file with the specified file path: {:s}",
66 path_to_tabular_file);
67 }
68
69 INFO("Found the tabular file: {:s}", path_to_tabular_file);
70
71 std::ifstream in(path_to_tabular_file);
72 if (!in)
73 {
74 OGS_FATAL("Couldn't open the tabular file: {:s}.",
75 path_to_tabular_file);
76 }
77
78 // read field names
79 std::string line;
80 std::getline(in, line);
81 std::vector<std::string> field_names;
82 boost::split(field_names, line, boost::is_any_of("\t "));
83
84 // categorize entry fields
85 std::vector<std::string> input_field_names;
86 std::copy_if(field_names.begin(), field_names.end(),
87 std::back_inserter(input_field_names),
88 [](auto const& field_name) -> bool
89 { return field_name.find("_new") == std::string::npos; });
90
91 // read table entries
92 std::map<std::string, std::vector<double>> tabular_data;
93 while (std::getline(in, line))
94 {
95 std::vector<std::string> field_data;
96 boost::split(field_data, line, boost::is_any_of("\t "));
97
98 assert(field_data.size() == field_names.size());
99 for (std::size_t field_id = 0; field_id < field_data.size(); ++field_id)
100 {
101 tabular_data[field_names[field_id]].push_back(
102 std::stod(field_data[field_id]));
103 }
104 }
105 in.close();
106
107 std::vector<Field> input_fields;
108 input_fields.reserve(input_field_names.size());
109 for (auto const& field_name : input_field_names)
110 {
111 auto pv = std::find_if(process_variables.begin(),
112 process_variables.end(),
113 [&field_name](auto const& v) -> bool
114 {
115 return v[0].get().getName() == field_name ||
116 v[0].get().getName() + "_prev" ==
117 field_name;
118 });
119
120 if (pv == process_variables.end())
121 {
122 OGS_FATAL(
123 "Not found field name {:s} in the group of process variables "
124 "defined in the project file.",
125 field_name);
126 }
127
128 auto const process_id =
129 static_cast<int>(std::distance(process_variables.cbegin(), pv));
130
131 auto seed_points = tabular_data[field_name];
132 BaseLib::makeVectorUnique(seed_points);
133
134 std::vector<std::vector<std::size_t>> point_id_groups;
135 for (auto const seed_point : seed_points)
136 {
137 auto const point_id_group =
138 getIndexVector(tabular_data[field_name], seed_point);
139 point_id_groups.push_back(point_id_group);
140 }
141
142 input_fields.emplace_back(point_id_groups, seed_points, field_name,
143 process_id);
144 }
145
146 return std::make_unique<LookupTable>(std::move(input_fields),
147 std::move(tabular_data));
148}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
bool IsFileExisting(const std::string &strFilename)
Returns true if given file exists.
Definition FileTools.cpp:23
std::string joinPaths(std::string const &pathA, std::string const &pathB)
void makeVectorUnique(std::vector< T > &v)
Definition Algorithm.h:173
std::vector< std::size_t > getIndexVector(std::vector< double > const &data, double const value)

References INFO(), BaseLib::IsFileExisting(), BaseLib::joinPaths(), BaseLib::makeVectorUnique(), and OGS_FATAL.

Referenced by createComponentTransportProcess().

◆ intersection()

void ProcessLib::ComponentTransport::intersection ( std::vector< std::size_t > & vec1,
std::vector< std::size_t > const & vec2 )
static

Definition at line 14 of file LookupTable.cpp.

16{
17 std::unordered_set<std::size_t> set(vec1.begin(), vec1.end());
18 vec1.clear();
19 for (auto const a : vec2)
20 {
21 if (set.contains(a))
22 {
23 vec1.push_back(a);
24 set.erase(a);
25 }
26 }
27}

Referenced by ProcessLib::ComponentTransport::LookupTable::getTableEntryID().