xref: /honee/src/velocity_gradient_projection.c (revision 0dee9b8e90ebabc518085db7cd2a8122f33364a4)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 /// @file
4 /// Functions for setting up and projecting the velocity gradient
5 
6 #include "../qfunctions/velocity_gradient_projection.h"
7 
8 #include <petscdmplex.h>
9 
10 #include <navierstokes.h>
11 
12 PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) {
13   PetscSection section;
14 
15   PetscFunctionBeginUser;
16   grad_velo_proj->num_comp = 9;  // 9 velocity gradient
17 
18   PetscCall(DMClone(user->dm, &grad_velo_proj->dm));
19   PetscCall(DMSetMatrixPreallocateSkip(grad_velo_proj->dm, PETSC_TRUE));
20   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
21 
22   PetscCall(
23       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));
24 
25   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
26   PetscCall(PetscSectionSetFieldName(section, 0, ""));
27   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
28   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
29   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
30   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
31   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
32   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
33   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
34   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
35   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
36   PetscFunctionReturn(PETSC_SUCCESS);
37 };
38 
39 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input,
40                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj) {
41   NodalProjectionData grad_velo_proj;
42   CeedOperator        op_rhs_assemble, op_mass;
43   CeedQFunction       qf_rhs_assemble, qf_mass;
44   CeedBasis           basis_grad_velo;
45   CeedElemRestriction elem_restr_grad_velo;
46   PetscInt            dim, label_value = 0, height = 0, dm_field = 0;
47   CeedInt             num_comp_x, num_comp_input, q_data_size;
48   DMLabel             domain_label = NULL;
49 
50   PetscFunctionBeginUser;
51   PetscCall(PetscNew(&grad_velo_proj));
52 
53   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
54 
55   // -- Get Pre-requisite things
56   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
57   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
58   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_input));
59   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
60   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &elem_restr_grad_velo));
61 
62   PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, domain_label, label_value, height, dm_field, &basis_grad_velo));
63 
64   // -- Build RHS operator
65   switch (state_var_input) {
66     case STATEVAR_PRIMITIVE:
67       PetscCallCeed(
68           ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble));
69       break;
70     case STATEVAR_CONSERVATIVE:
71       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc,
72                                                       &qf_rhs_assemble));
73       break;
74     case STATEVAR_ENTROPY:
75       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Entropy, VelocityGradientProjectionRHS_Entropy_loc,
76                                                       &qf_rhs_assemble));
77       break;
78   }
79 
80   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfctx));
81   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_input, CEED_EVAL_INTERP));
82   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_input * dim, CEED_EVAL_GRAD));
83   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE));
84   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP));
85   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP));
86 
87   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble));
88   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE));
89   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", elem_restr_input, basis_input, CEED_VECTOR_ACTIVE));
90   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
91   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
92   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
93 
94   PetscCall(OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, NULL, NULL, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
95 
96   // -- Build Mass operator
97   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
98   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
99   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
100   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
101   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE));
102 
103   {  // -- Setup KSP for L^2 projection with lumped mass operator
104     Mat      mat_mass;
105     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
106 
107     PetscCall(MatCreateCeed(grad_velo_proj->dm, grad_velo_proj->dm, op_mass, NULL, &mat_mass));
108 
109     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
110     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
111     {
112       PC pc;
113       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
114       PetscCall(PCSetType(pc, PCJACOBI));
115       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
116       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
117     }
118     PetscCall(KSPSetFromOptions_WithMatCeed(grad_velo_proj->ksp, mat_mass));
119     PetscCall(MatDestroy(&mat_mass));
120   }
121 
122   *pgrad_velo_proj = grad_velo_proj;
123 
124   PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo));
125   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
126   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble));
127   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
128   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble));
129   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
130   PetscFunctionReturn(PETSC_SUCCESS);
131 }
132 
133 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient) {
134   OperatorApplyContext l2_rhs_ctx = grad_velo_proj->l2_rhs_ctx;
135 
136   PetscFunctionBeginUser;
137   PetscCall(PetscLogEventBegin(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0));
138   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
139 
140   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
141   PetscCall(PetscLogEventEnd(FLUIDS_VelocityGradientProjection, Q_loc, VelocityGradient, 0, 0));
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144