xref: /honee/src/velocity_gradient_projection.c (revision 58e1cbfda2c067b7e8bf8571e993e896d411c8f5)
1457a5831SJames Wright // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors.
2457a5831SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3457a5831SJames Wright //
4457a5831SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5457a5831SJames Wright //
6457a5831SJames Wright // This file is part of CEED:  http://github.com/ceed
7457a5831SJames Wright /// @file
8457a5831SJames Wright /// Functions for setting up and projecting the velocity gradient
9457a5831SJames Wright 
10457a5831SJames Wright #include "../qfunctions/velocity_gradient_projection.h"
11457a5831SJames Wright 
12457a5831SJames Wright #include <petscdmplex.h>
13457a5831SJames Wright 
14457a5831SJames Wright #include "../navierstokes.h"
15457a5831SJames Wright 
16457a5831SJames Wright PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) {
17457a5831SJames Wright   PetscSection section;
18457a5831SJames Wright 
19457a5831SJames Wright   PetscFunctionBeginUser;
20457a5831SJames Wright   grad_velo_proj->num_comp = 9;  // 9 velocity gradient
21457a5831SJames Wright 
22457a5831SJames Wright   PetscCall(DMClone(user->dm, &grad_velo_proj->dm));
23457a5831SJames Wright   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
24457a5831SJames Wright 
25da4ca0cfSJames Wright   PetscCall(
26da4ca0cfSJames 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));
27457a5831SJames Wright 
28457a5831SJames Wright   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
29457a5831SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
30457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
31457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
32457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
33457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
34457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
35457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
36457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
37457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
38457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
39d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
40457a5831SJames Wright };
419ab09d51SJames Wright 
429ab09d51SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
439ab09d51SJames Wright   OperatorApplyContext mass_matop_ctx;
449ab09d51SJames Wright   CeedOperator         op_rhs_assemble, op_mass;
459ab09d51SJames Wright   CeedQFunction        qf_rhs_assemble, qf_mass;
469ab09d51SJames Wright   CeedBasis            basis_grad_velo;
479ab09d51SJames Wright   CeedElemRestriction  elem_restr_grad_velo;
4815c18037SJames Wright   PetscInt             dim, label_value = 0, height = 0, dm_field = 0;
4967263decSKenneth E. Jansen   CeedInt              num_comp_x, num_comp_q, q_data_size;
5015c18037SJames Wright   DMLabel              domain_label = NULL;
519ab09d51SJames Wright 
529ab09d51SJames Wright   PetscFunctionBeginUser;
539ab09d51SJames Wright   PetscCall(PetscNew(&user->grad_velo_proj));
549ab09d51SJames Wright   NodalProjectionData grad_velo_proj = user->grad_velo_proj;
559ab09d51SJames Wright 
569ab09d51SJames Wright   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
579ab09d51SJames Wright 
589ab09d51SJames Wright   // -- Get Pre-requisite things
599ab09d51SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
60b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
61b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
62b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
6315c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &elem_restr_grad_velo));
649ab09d51SJames Wright 
6515c18037SJames Wright   PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &basis_grad_velo));
669ab09d51SJames Wright 
679ab09d51SJames Wright   // -- Build RHS operator
689ab09d51SJames Wright   switch (user->phys->state_var) {
699ab09d51SJames Wright     case STATEVAR_PRIMITIVE:
70b4c37c5cSJames Wright       PetscCallCeed(
71b4c37c5cSJames Wright           ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble));
729ab09d51SJames Wright       break;
739ab09d51SJames Wright     case STATEVAR_CONSERVATIVE:
74b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc,
75b4c37c5cSJames Wright                                                       &qf_rhs_assemble));
769ab09d51SJames Wright       break;
779ab09d51SJames Wright     default:
789ab09d51SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable");
799ab09d51SJames Wright   }
809ab09d51SJames Wright 
81b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context));
82b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_q, CEED_EVAL_INTERP));
83b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
84b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE));
85b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP));
86b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP));
879ab09d51SJames Wright 
88b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble));
89b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
90b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
91*58e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
92b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
93b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
949ab09d51SJames Wright 
9549f5db74SJames Wright   PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
969ab09d51SJames Wright 
979ab09d51SJames Wright   // -- Build Mass operator
989ab09d51SJames Wright   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
99b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
100b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
101*58e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
102b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
1039ab09d51SJames Wright 
1049ab09d51SJames Wright   {  // -- Setup KSP for L^2 projection with lumped mass operator
1059ab09d51SJames Wright     Mat      mat_mass;
1069ab09d51SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
1079ab09d51SJames Wright 
10849f5db74SJames Wright     PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, grad_velo_proj->dm, ceed, op_mass, NULL, NULL, NULL, NULL, &mass_matop_ctx));
10949f5db74SJames Wright     PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass));
1109ab09d51SJames Wright 
1119ab09d51SJames Wright     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
1129ab09d51SJames Wright     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
1139ab09d51SJames Wright     {
1149ab09d51SJames Wright       PC pc;
1159ab09d51SJames Wright       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
1169ab09d51SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
1179ab09d51SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1189ab09d51SJames Wright       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
1199ab09d51SJames Wright       // TODO Not sure if the option below are necessary
1209ab09d51SJames Wright       PetscCall(KSPSetConvergenceTest(grad_velo_proj->ksp, KSPConvergedSkip, NULL, NULL));
1219ab09d51SJames Wright     }
1229ab09d51SJames Wright     PetscCall(KSPSetOperators(grad_velo_proj->ksp, mat_mass, mat_mass));
1239ab09d51SJames Wright     PetscCall(KSPSetFromOptions(grad_velo_proj->ksp));
1249ab09d51SJames Wright   }
1259ab09d51SJames Wright 
126b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo));
127b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
128b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble));
129b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
130b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble));
131b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
132d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1339ab09d51SJames Wright }
134a2401be5SJames Wright 
135a2401be5SJames Wright PetscErrorCode VelocityGradientProjectionApply(User user, Vec Q_loc, Vec VelocityGradient) {
136a2401be5SJames Wright   NodalProjectionData  grad_velo_proj = user->grad_velo_proj;
137a2401be5SJames Wright   OperatorApplyContext l2_rhs_ctx     = grad_velo_proj->l2_rhs_ctx;
138a2401be5SJames Wright 
139a2401be5SJames Wright   PetscFunctionBeginUser;
140a2401be5SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
141a2401be5SJames Wright 
142a2401be5SJames Wright   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
143d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
144a2401be5SJames Wright }
145