xref: /libCEED/examples/fluids/src/velocity_gradient_projection.c (revision 3737f8325a1e09ea5db4e0ce05b04230e5d3648e)
1999ff5c7SJames Wright // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors.
2999ff5c7SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3999ff5c7SJames Wright //
4999ff5c7SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5999ff5c7SJames Wright //
6999ff5c7SJames Wright // This file is part of CEED:  http://github.com/ceed
7999ff5c7SJames Wright /// @file
8999ff5c7SJames Wright /// Functions for setting up and projecting the velocity gradient
9999ff5c7SJames Wright 
10999ff5c7SJames Wright #include "../qfunctions/velocity_gradient_projection.h"
11999ff5c7SJames Wright 
12999ff5c7SJames Wright #include <petscdmplex.h>
13999ff5c7SJames Wright 
14999ff5c7SJames Wright #include "../navierstokes.h"
15999ff5c7SJames Wright 
16999ff5c7SJames Wright PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) {
17999ff5c7SJames Wright   PetscSection section;
18999ff5c7SJames Wright 
19999ff5c7SJames Wright   PetscFunctionBeginUser;
20999ff5c7SJames Wright   grad_velo_proj->num_comp = 9;  // 9 velocity gradient
21999ff5c7SJames Wright 
22999ff5c7SJames Wright   PetscCall(DMClone(user->dm, &grad_velo_proj->dm));
23999ff5c7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
24999ff5c7SJames Wright 
25d68a66c4SJames Wright   PetscCall(
26d68a66c4SJames Wright       DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &grad_velo_proj->num_comp, grad_velo_proj->dm));
27999ff5c7SJames Wright 
28999ff5c7SJames Wright   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
29999ff5c7SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
30999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
31999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
32999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
33999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
34999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
35999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
36999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
37999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
38999ff5c7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
39ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
40999ff5c7SJames Wright };
413909d146SJames Wright 
427618d7b7SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, StateVariable state_var_input,
437618d7b7SJames Wright                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj) {
447618d7b7SJames Wright   NodalProjectionData  grad_velo_proj;
453909d146SJames Wright   OperatorApplyContext mass_matop_ctx;
463909d146SJames Wright   CeedOperator         op_rhs_assemble, op_mass;
473909d146SJames Wright   CeedQFunction        qf_rhs_assemble, qf_mass;
483909d146SJames Wright   CeedBasis            basis_grad_velo;
493909d146SJames Wright   CeedElemRestriction  elem_restr_grad_velo;
50bb85d312SJames Wright   PetscInt             dim, label_value = 0, height = 0, dm_field = 0;
517618d7b7SJames Wright   CeedInt              num_comp_x, num_comp_input, q_data_size;
52bb85d312SJames Wright   DMLabel              domain_label = NULL;
533909d146SJames Wright 
543909d146SJames Wright   PetscFunctionBeginUser;
557618d7b7SJames Wright   PetscCall(PetscNew(&grad_velo_proj));
563909d146SJames Wright 
573909d146SJames Wright   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
583909d146SJames Wright 
593909d146SJames Wright   // -- Get Pre-requisite things
603909d146SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
61a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
627618d7b7SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_input));
63a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
64bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &elem_restr_grad_velo));
653909d146SJames Wright 
66bb85d312SJames Wright   PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &basis_grad_velo));
673909d146SJames Wright 
683909d146SJames Wright   // -- Build RHS operator
697618d7b7SJames Wright   switch (state_var_input) {
703909d146SJames Wright     case STATEVAR_PRIMITIVE:
71a424bcd0SJames Wright       PetscCallCeed(
72a424bcd0SJames Wright           ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble));
733909d146SJames Wright       break;
743909d146SJames Wright     case STATEVAR_CONSERVATIVE:
75a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc,
76a424bcd0SJames Wright                                                       &qf_rhs_assemble));
773909d146SJames Wright       break;
783909d146SJames Wright     default:
793909d146SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable");
803909d146SJames Wright   }
813909d146SJames Wright 
82a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context));
837618d7b7SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_input, CEED_EVAL_INTERP));
847618d7b7SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_input * dim, CEED_EVAL_GRAD));
85a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE));
86a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP));
87a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP));
883909d146SJames Wright 
89a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble));
907618d7b7SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE));
917618d7b7SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE));
92356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
93a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
94a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
953909d146SJames Wright 
96e66242d1SJames Wright   PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
973909d146SJames Wright 
983909d146SJames Wright   // -- Build Mass operator
993909d146SJames Wright   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
100a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
101a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
102356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
103a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
1043909d146SJames Wright 
1053909d146SJames Wright   {  // -- Setup KSP for L^2 projection with lumped mass operator
1063909d146SJames Wright     Mat      mat_mass;
1073909d146SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
1083909d146SJames Wright 
109e66242d1SJames Wright     PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, grad_velo_proj->dm, ceed, op_mass, NULL, NULL, NULL, NULL, &mass_matop_ctx));
110e66242d1SJames Wright     PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass));
1113909d146SJames Wright 
1123909d146SJames Wright     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
1133909d146SJames Wright     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
1143909d146SJames Wright     {
1153909d146SJames Wright       PC pc;
1163909d146SJames Wright       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
1173909d146SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
1183909d146SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1193909d146SJames Wright       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
1203909d146SJames Wright       // TODO Not sure if the option below are necessary
1213909d146SJames Wright       PetscCall(KSPSetConvergenceTest(grad_velo_proj->ksp, KSPConvergedSkip, NULL, NULL));
1223909d146SJames Wright     }
1233909d146SJames Wright     PetscCall(KSPSetOperators(grad_velo_proj->ksp, mat_mass, mat_mass));
1243909d146SJames Wright     PetscCall(KSPSetFromOptions(grad_velo_proj->ksp));
1253909d146SJames Wright   }
1263909d146SJames Wright 
1277618d7b7SJames Wright   *pgrad_velo_proj = grad_velo_proj;
1287618d7b7SJames Wright 
129a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo));
130a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
131a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble));
132a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
133a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble));
134a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
135ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1363909d146SJames Wright }
137f84dcd56SJames Wright 
1387618d7b7SJames Wright PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient) {
139f84dcd56SJames Wright   OperatorApplyContext l2_rhs_ctx = grad_velo_proj->l2_rhs_ctx;
140f84dcd56SJames Wright 
141f84dcd56SJames Wright   PetscFunctionBeginUser;
142*3737f832SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0));
143f84dcd56SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
144f84dcd56SJames Wright 
145f84dcd56SJames Wright   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
146*3737f832SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0));
147ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
148f84dcd56SJames Wright }
149