xref: /honee/src/velocity_gradient_projection.c (revision 9ab09d51ac9d937f9ba51ebe6202d581c3c26d67)
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 
28457a5831SJames Wright   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, grad_velo_proj->num_comp, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
29457a5831SJames Wright   PetscCall(PetscObjectSetName((PetscObject)fe, "Velocity Gradient Projection"));
30457a5831SJames Wright   PetscCall(DMAddField(grad_velo_proj->dm, NULL, (PetscObject)fe));
31457a5831SJames Wright   PetscCall(DMCreateDS(grad_velo_proj->dm));
32457a5831SJames Wright   PetscCall(DMPlexSetClosurePermutationTensor(grad_velo_proj->dm, PETSC_DETERMINE, NULL));
33457a5831SJames Wright 
34457a5831SJames Wright   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
35457a5831SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
36457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
37457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
38457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
39457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
40457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
41457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
42457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
43457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
44457a5831SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
45457a5831SJames Wright 
46457a5831SJames Wright   PetscCall(PetscFEDestroy(&fe));
47457a5831SJames Wright   PetscFunctionReturn(0);
48457a5831SJames Wright };
49*9ab09d51SJames Wright 
50*9ab09d51SJames Wright PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
51*9ab09d51SJames Wright   OperatorApplyContext mass_matop_ctx;
52*9ab09d51SJames Wright   CeedOperator         op_rhs_assemble, op_mass;
53*9ab09d51SJames Wright   CeedQFunction        qf_rhs_assemble, qf_mass;
54*9ab09d51SJames Wright   CeedBasis            basis_grad_velo;
55*9ab09d51SJames Wright   CeedVector           q_ceed, rhs_ceed, mass_output;
56*9ab09d51SJames Wright   CeedElemRestriction  elem_restr_grad_velo;
57*9ab09d51SJames Wright   PetscInt             dim, num_comp_x, num_comp_q, q_data_size, num_qpts_1d, num_nodes_1d;
58*9ab09d51SJames Wright 
59*9ab09d51SJames Wright   PetscFunctionBeginUser;
60*9ab09d51SJames Wright   PetscCall(PetscNew(&user->grad_velo_proj));
61*9ab09d51SJames Wright   NodalProjectionData grad_velo_proj = user->grad_velo_proj;
62*9ab09d51SJames Wright 
63*9ab09d51SJames Wright   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
64*9ab09d51SJames Wright 
65*9ab09d51SJames Wright   // -- Get Pre-requisite things
66*9ab09d51SJames Wright   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
67*9ab09d51SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
68*9ab09d51SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
69*9ab09d51SJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &num_qpts_1d);
70*9ab09d51SJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_q, &num_nodes_1d);
71*9ab09d51SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
72*9ab09d51SJames Wright 
73*9ab09d51SJames Wright   PetscCall(GetRestrictionForDomain(ceed, grad_velo_proj->dm, 0, 0, 0, num_qpts_1d, q_data_size, &elem_restr_grad_velo, NULL, NULL));
74*9ab09d51SJames Wright 
75*9ab09d51SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, grad_velo_proj->num_comp, num_nodes_1d, num_qpts_1d, CEED_GAUSS, &basis_grad_velo);
76*9ab09d51SJames Wright 
77*9ab09d51SJames Wright   // -- Build RHS operator
78*9ab09d51SJames Wright   switch (user->phys->state_var) {
79*9ab09d51SJames Wright     case STATEVAR_PRIMITIVE:
80*9ab09d51SJames Wright       CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble);
81*9ab09d51SJames Wright       break;
82*9ab09d51SJames Wright     case STATEVAR_CONSERVATIVE:
83*9ab09d51SJames Wright       CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc, &qf_rhs_assemble);
84*9ab09d51SJames Wright       break;
85*9ab09d51SJames Wright     default:
86*9ab09d51SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable");
87*9ab09d51SJames Wright   }
88*9ab09d51SJames Wright 
89*9ab09d51SJames Wright   CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context);
90*9ab09d51SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_q, CEED_EVAL_INTERP);
91*9ab09d51SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
92*9ab09d51SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "qdata", q_data_size, CEED_EVAL_NONE);
93*9ab09d51SJames Wright   CeedQFunctionAddInput(qf_rhs_assemble, "x", num_comp_x, CEED_EVAL_INTERP);
94*9ab09d51SJames Wright   CeedQFunctionAddOutput(qf_rhs_assemble, "velocity gradient", grad_velo_proj->num_comp, CEED_EVAL_INTERP);
95*9ab09d51SJames Wright 
96*9ab09d51SJames Wright   CeedOperatorCreate(ceed, qf_rhs_assemble, NULL, NULL, &op_rhs_assemble);
97*9ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
98*9ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
99*9ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
100*9ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
101*9ab09d51SJames Wright   CeedOperatorSetField(op_rhs_assemble, "velocity gradient", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
102*9ab09d51SJames Wright 
103*9ab09d51SJames Wright   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q_ceed, NULL);
104*9ab09d51SJames Wright   CeedElemRestrictionCreateVector(elem_restr_grad_velo, &rhs_ceed, NULL);
105*9ab09d51SJames Wright 
106*9ab09d51SJames Wright   PetscCall(
107*9ab09d51SJames Wright       OperatorApplyContextCreate(user->dm, grad_velo_proj->dm, ceed, op_rhs_assemble, q_ceed, rhs_ceed, NULL, NULL, &grad_velo_proj->l2_rhs_ctx));
108*9ab09d51SJames Wright 
109*9ab09d51SJames Wright   // -- Build Mass operator
110*9ab09d51SJames Wright   PetscCall(CreateMassQFunction(ceed, grad_velo_proj->num_comp, q_data_size, &qf_mass));
111*9ab09d51SJames Wright   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
112*9ab09d51SJames Wright   CeedOperatorSetField(op_mass, "u", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
113*9ab09d51SJames Wright   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
114*9ab09d51SJames Wright   CeedOperatorSetField(op_mass, "v", elem_restr_grad_velo, basis_grad_velo, CEED_VECTOR_ACTIVE);
115*9ab09d51SJames Wright 
116*9ab09d51SJames Wright   {  // -- Setup KSP for L^2 projection with lumped mass operator
117*9ab09d51SJames Wright     PetscInt l_size, g_size;
118*9ab09d51SJames Wright     Mat      mat_mass;
119*9ab09d51SJames Wright     VecType  vec_type;
120*9ab09d51SJames Wright     Vec      M_inv;
121*9ab09d51SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)grad_velo_proj->dm);
122*9ab09d51SJames Wright 
123*9ab09d51SJames Wright     PetscCall(DMGetGlobalVector(grad_velo_proj->dm, &M_inv));
124*9ab09d51SJames Wright     PetscCall(VecGetLocalSize(M_inv, &l_size));
125*9ab09d51SJames Wright     PetscCall(VecGetSize(M_inv, &g_size));
126*9ab09d51SJames Wright     PetscCall(VecGetType(M_inv, &vec_type));
127*9ab09d51SJames Wright     PetscCall(DMRestoreGlobalVector(grad_velo_proj->dm, &M_inv));
128*9ab09d51SJames Wright 
129*9ab09d51SJames Wright     CeedElemRestrictionCreateVector(elem_restr_grad_velo, &mass_output, NULL);
130*9ab09d51SJames Wright     PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, grad_velo_proj->dm, ceed, op_mass, rhs_ceed, mass_output, NULL, NULL, &mass_matop_ctx));
131*9ab09d51SJames Wright     CeedVectorDestroy(&mass_output);
132*9ab09d51SJames Wright 
133*9ab09d51SJames Wright     PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, mass_matop_ctx, &mat_mass));
134*9ab09d51SJames Wright     PetscCall(MatShellSetContextDestroy(mat_mass, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
135*9ab09d51SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed));
136*9ab09d51SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
137*9ab09d51SJames Wright     PetscCall(MatShellSetVecType(mat_mass, vec_type));
138*9ab09d51SJames Wright 
139*9ab09d51SJames Wright     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
140*9ab09d51SJames Wright     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
141*9ab09d51SJames Wright     {
142*9ab09d51SJames Wright       PC pc;
143*9ab09d51SJames Wright       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
144*9ab09d51SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
145*9ab09d51SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
146*9ab09d51SJames Wright       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
147*9ab09d51SJames Wright       // TODO Not sure if the option below are necessary
148*9ab09d51SJames Wright       PetscCall(KSPSetConvergenceTest(grad_velo_proj->ksp, KSPConvergedSkip, NULL, NULL));
149*9ab09d51SJames Wright     }
150*9ab09d51SJames Wright     PetscCall(KSPSetOperators(grad_velo_proj->ksp, mat_mass, mat_mass));
151*9ab09d51SJames Wright     PetscCall(KSPSetFromOptions(grad_velo_proj->ksp));
152*9ab09d51SJames Wright   }
153*9ab09d51SJames Wright 
154*9ab09d51SJames Wright   CeedVectorDestroy(&q_ceed);
155*9ab09d51SJames Wright   CeedVectorDestroy(&rhs_ceed);
156*9ab09d51SJames Wright   CeedBasisDestroy(&basis_grad_velo);
157*9ab09d51SJames Wright   CeedElemRestrictionDestroy(&elem_restr_grad_velo);
158*9ab09d51SJames Wright   CeedQFunctionDestroy(&qf_rhs_assemble);
159*9ab09d51SJames Wright   CeedQFunctionDestroy(&qf_mass);
160*9ab09d51SJames Wright   CeedOperatorDestroy(&op_rhs_assemble);
161*9ab09d51SJames Wright   CeedOperatorDestroy(&op_mass);
162*9ab09d51SJames Wright   PetscFunctionReturn(0);
163*9ab09d51SJames Wright }
164