Skip to content

Commit

Permalink
account for the fluid velocity for outflow boundaries
Browse files Browse the repository at this point in the history
  • Loading branch information
danieldouglas92 committed Oct 14, 2024
1 parent 18b35da commit b52900a
Showing 1 changed file with 17 additions and 6 deletions.
23 changes: 17 additions & 6 deletions source/simulator/helper_functions.cc
Original file line number Diff line number Diff line change
Expand Up @@ -2291,12 +2291,12 @@ namespace aspect
update_quadrature_points | update_JxW_values);

std::vector<Tensor<1,dim>> face_current_velocity_values (fe_face_values.n_quadrature_points);
std::vector<Tensor<1,dim>> face_current_fluid_velocity_values(fe_face_values.n_quadrature_points);
std::vector<Tensor<1,dim>> face_current_fluid_velocity_values (fe_face_values.n_quadrature_points);

MaterialModel::MaterialModelInputs<dim> in(fe_face_values.n_quadrature_points, introspection.n_compositional_fields);
MaterialModel::MaterialModelOutputs<dim> out(fe_face_values.n_quadrature_points, introspection.n_compositional_fields);
MeltHandler<dim>::create_material_model_outputs(out);
MaterialModel::MeltOutputs<dim> *fluid_out = nullptr;
MaterialModel::MeltOutputs<dim> *fluid_out = out.template get_additional_output<MaterialModel::MeltOutputs<dim>>();;

const auto &tangential_velocity_boundaries =
boundary_velocity_manager.get_tangential_boundary_velocity_indicators();
Expand Down Expand Up @@ -2329,8 +2329,11 @@ namespace aspect

// ... check if the face is an outflow boundary by integrating the normal velocities
// (flux through the boundary) as: int u*n ds = Sum_q u(x_q)*n(x_q) JxW(x_q)...
double integrated_flow = 0;
// do this for the solid velocity, the darcy velocity, and/or the fluid velocity
// if applicable.
double integrated_solid_flow = 0;
double integrated_darcy_flow = 0;
double integrated_fluid_flow = 0;

for (unsigned int q=0; q<fe_face_values.n_quadrature_points; ++q)
{
Expand All @@ -2356,12 +2359,20 @@ namespace aspect
fe_face_values.JxW(q);
}

integrated_flow += (boundary_velocity * fe_face_values.normal_vector(q)) *
fe_face_values.JxW(q);
if (parameters.include_melt_transport == true)
{
const FEValuesExtractors::Vector ex_u_f = introspection.variable("fluid velocity").extractor_vector();
fe_face_values[ex_u_f].get_function_values(current_linearization_point, face_current_fluid_velocity_values);
integrated_fluid_flow += (face_current_fluid_velocity_values[q] * fe_face_values.normal_vector(q)) *
fe_face_values.JxW(q);
}

integrated_solid_flow += (boundary_velocity * fe_face_values.normal_vector(q)) *
fe_face_values.JxW(q);
}

// ... and change the boundary id of any outflow boundary faces.
if (integrated_flow > 0 || integrated_darcy_flow > 0)
if (integrated_solid_flow > 0 || integrated_darcy_flow > 0 || integrated_fluid_flow > 0)
face->set_boundary_id(face->boundary_id() + offset);
}
}
Expand Down

0 comments on commit b52900a

Please sign in to comment.