xref: /honee/src/velocity_gradient_projection.c (revision b4c37c5c26bd208d7f474cbd9aa0850e82e1a854)
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   PetscFE      fe;
18457a5831SJames Wright   PetscSection section;
19457a5831SJames Wright   PetscInt     dim;
20457a5831SJames Wright 
21457a5831SJames Wright   PetscFunctionBeginUser;
22457a5831SJames Wright   grad_velo_proj->num_comp = 9;  // 9 velocity gradient
23457a5831SJames Wright 
24457a5831SJames Wright   PetscCall(DMClone(user->dm, &grad_velo_proj->dm));
25457a5831SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
26457a5831SJames Wright   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
27457a5831SJames Wright 
2867263decSKenneth E. Jansen   PetscInt q_order = user->app_ctx->degree + user->app_ctx->q_extra;
2967263decSKenneth E. Jansen   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, grad_velo_proj->num_comp, PETSC_FALSE, degree, q_order, &fe));
30457a5831SJames Wright   PetscCall(PetscObjectSetName((PetscObject)fe, "Velocity Gradient Projection"));
31457a5831SJames Wright   PetscCall(DMAddField(grad_velo_proj->dm, NULL, (PetscObject)fe));
32457a5831SJames Wright   PetscCall(DMCreateDS(grad_velo_proj->dm));
33457a5831SJames Wright   PetscCall(DMPlexSetClosurePermutationTensor(grad_velo_proj->dm, PETSC_DETERMINE, NULL));
34457a5831SJames Wright 
35457a5831SJames Wright   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
36457a5831SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
37457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
38457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
39457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
40457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
41457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
42457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
43457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
44457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
45457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
46457a5831SJames Wright 
47457a5831SJames Wright   PetscCall(PetscFEDestroy(&fe));
48d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
49457a5831SJames Wright };
509ab09d51SJames Wright 
519ab09d51SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
529ab09d51SJames Wright   OperatorApplyContext mass_matop_ctx;
539ab09d51SJames Wright   CeedOperator         op_rhs_assemble, op_mass;
549ab09d51SJames Wright   CeedQFunction        qf_rhs_assemble, qf_mass;
559ab09d51SJames Wright   CeedBasis            basis_grad_velo;
569ab09d51SJames Wright   CeedElemRestriction  elem_restr_grad_velo;
57defe8520SJames Wright   PetscInt             dim;
5867263decSKenneth E. Jansen   CeedInt              num_comp_x, num_comp_q, q_data_size;
599ab09d51SJames Wright 
609ab09d51SJames Wright   PetscFunctionBeginUser;
619ab09d51SJames Wright   PetscCall(PetscNew(&user->grad_velo_proj));
629ab09d51SJames Wright   NodalProjectionData grad_velo_proj = user->grad_velo_proj;
639ab09d51SJames Wright 
649ab09d51SJames Wright   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
659ab09d51SJames Wright 
669ab09d51SJames Wright   // -- Get Pre-requisite things
679ab09d51SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
68*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
69*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
70*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
7167263decSKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, grad_velo_proj->dm, 0, 0, 0, 0, -1, 0, &elem_restr_grad_velo, NULL, NULL));
729ab09d51SJames Wright 
7367263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, 0, 0, 0, 0, &basis_grad_velo));
749ab09d51SJames Wright 
759ab09d51SJames Wright   // -- Build RHS operator
769ab09d51SJames Wright   switch (user->phys->state_var) {
779ab09d51SJames Wright     case STATEVAR_PRIMITIVE:
78*b4c37c5cSJames Wright       PetscCallCeed(
79*b4c37c5cSJames Wright           ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble));
809ab09d51SJames Wright       break;
819ab09d51SJames Wright     case STATEVAR_CONSERVATIVE:
82*b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc,
83*b4c37c5cSJames Wright                                                       &qf_rhs_assemble));
849ab09d51SJames Wright       break;
859ab09d51SJames Wright     default:
869ab09d51SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable");
879ab09d51SJames Wright   }
889ab09d51SJames Wright 
89*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context));
90*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_q, CEED_EVAL_INTERP));
91*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
92*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE));
93*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP));
94*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP));
959ab09d51SJames Wright 
96*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble));
97*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
98*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
99*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
100*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
101*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
1029ab09d51SJames Wright 
10349f5db74SJames Wright   PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
1049ab09d51SJames Wright 
1059ab09d51SJames Wright   // -- Build Mass operator
1069ab09d51SJames Wright   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
107*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
108*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
109*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
110*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
1119ab09d51SJames Wright 
1129ab09d51SJames Wright   {  // -- Setup KSP for L^2 projection with lumped mass operator
1139ab09d51SJames Wright     Mat      mat_mass;
1149ab09d51SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
1159ab09d51SJames Wright 
11649f5db74SJames Wright     PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, grad_velo_proj->dm, ceed, op_mass, NULL, NULL, NULL, NULL, &mass_matop_ctx));
11749f5db74SJames Wright     PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass));
1189ab09d51SJames Wright 
1199ab09d51SJames Wright     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
1209ab09d51SJames Wright     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
1219ab09d51SJames Wright     {
1229ab09d51SJames Wright       PC pc;
1239ab09d51SJames Wright       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
1249ab09d51SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
1259ab09d51SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1269ab09d51SJames Wright       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
1279ab09d51SJames Wright       // TODO Not sure if the option below are necessary
1289ab09d51SJames Wright       PetscCall(KSPSetConvergenceTest(grad_velo_proj->ksp, KSPConvergedSkip, NULL, NULL));
1299ab09d51SJames Wright     }
1309ab09d51SJames Wright     PetscCall(KSPSetOperators(grad_velo_proj->ksp, mat_mass, mat_mass));
1319ab09d51SJames Wright     PetscCall(KSPSetFromOptions(grad_velo_proj->ksp));
1329ab09d51SJames Wright   }
1339ab09d51SJames Wright 
134*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo));
135*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
136*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble));
137*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
138*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble));
139*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
140d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1419ab09d51SJames Wright }
142a2401be5SJames Wright 
143a2401be5SJames Wright PetscErrorCode VelocityGradientProjectionApply(User user, Vec Q_loc, Vec VelocityGradient) {
144a2401be5SJames Wright   NodalProjectionData  grad_velo_proj = user->grad_velo_proj;
145a2401be5SJames Wright   OperatorApplyContext l2_rhs_ctx     = grad_velo_proj->l2_rhs_ctx;
146a2401be5SJames Wright 
147a2401be5SJames Wright   PetscFunctionBeginUser;
148a2401be5SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
149a2401be5SJames Wright 
150a2401be5SJames Wright   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
151a2401be5SJames Wright 
152d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
153a2401be5SJames Wright }
154