xref: /libCEED/examples/fluids/src/velocity_gradient_projection.c (revision 524ffdfd70d97d4e960fc9822927ee4b8f8d8e10)
1 // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 /// @file
8 /// Functions for setting up and projecting the velocity gradient
9 
10 #include "../qfunctions/velocity_gradient_projection.h"
11 
12 #include <petscdmplex.h>
13 
14 #include "../navierstokes.h"
15 
16 PetscErrorCode VelocityGradientProjectionCreateDM(NodalProjectionData grad_velo_proj, User user, PetscInt degree) {
17   PetscSection section;
18 
19   PetscFunctionBeginUser;
20   grad_velo_proj->num_comp = 9;  // 9 velocity gradient
21 
22   PetscCall(DMClone(user->dm, &grad_velo_proj->dm));
23   PetscCall(PetscObjectSetName((PetscObject)grad_velo_proj->dm, "Velocity Gradient Projection"));
24 
25   PetscCall(
26       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));
27 
28   PetscCall(DMGetLocalSection(grad_velo_proj->dm, &section));
29   PetscCall(PetscSectionSetFieldName(section, 0, ""));
30   PetscCall(PetscSectionSetComponentName(section, 0, 0, "VelocityGradientXX"));
31   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityGradientXY"));
32   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityGradientXZ"));
33   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityGradientYX"));
34   PetscCall(PetscSectionSetComponentName(section, 0, 4, "VelocityGradientYY"));
35   PetscCall(PetscSectionSetComponentName(section, 0, 5, "VelocityGradientYZ"));
36   PetscCall(PetscSectionSetComponentName(section, 0, 6, "VelocityGradientZX"));
37   PetscCall(PetscSectionSetComponentName(section, 0, 7, "VelocityGradientZY"));
38   PetscCall(PetscSectionSetComponentName(section, 0, 8, "VelocityGradientZZ"));
39   PetscFunctionReturn(PETSC_SUCCESS);
40 };
41 
42 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
43   OperatorApplyContext mass_matop_ctx;
44   CeedOperator         op_rhs_assemble, op_mass;
45   CeedQFunction        qf_rhs_assemble, qf_mass;
46   CeedBasis            basis_grad_velo;
47   CeedElemRestriction  elem_restr_grad_velo;
48   PetscInt             dim;
49   CeedInt              num_comp_x, num_comp_q, q_data_size;
50 
51   PetscFunctionBeginUser;
52   PetscCall(PetscNew(&user->grad_velo_proj));
53   NodalProjectionData grad_velo_proj = user->grad_velo_proj;
54 
55   PetscCall(VelocityGradientProjectionCreateDM(grad_velo_proj, user, user->app_ctx->degree));
56 
57   // -- Get Pre-requisite things
58   PetscCall(DMGetDimension(grad_velo_proj->dm, &dim));
59   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
60   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
61   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
62   PetscCall(GetRestrictionForDomain(ceed, grad_velo_proj->dm, 0, 0, 0, 0, -1, 0, &elem_restr_grad_velo, NULL, NULL));
63 
64   PetscCall(CreateBasisFromPlex(ceed, grad_velo_proj->dm, 0, 0, 0, 0, &basis_grad_velo));
65 
66   // -- Build RHS operator
67   switch (user->phys->state_var) {
68     case STATEVAR_PRIMITIVE:
69       PetscCallCeed(
70           ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Prim, VelocityGradientProjectionRHS_Prim_loc, &qf_rhs_assemble));
71       break;
72     case STATEVAR_CONSERVATIVE:
73       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, VelocityGradientProjectionRHS_Conserv, VelocityGradientProjectionRHS_Conserv_loc,
74                                                       &qf_rhs_assemble));
75       break;
76     default:
77       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No velocity gradient projection QFunction for chosen state variable");
78   }
79 
80   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_assemble, problem->apply_vol_ifunction.qfunction_context));
81   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "q", num_comp_q, CEED_EVAL_INTERP));
82   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_assemble, "Grad_q", num_comp_q * 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", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
89   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
90   PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_assemble, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, 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_COLLOCATED, 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(OperatorApplyContextCreate(grad_velo_proj->dm, grad_velo_proj->dm, ceed, op_mass, NULL, NULL, NULL, NULL, &mass_matop_ctx));
108     PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass));
109 
110     PetscCall(KSPCreate(comm, &grad_velo_proj->ksp));
111     PetscCall(KSPSetOptionsPrefix(grad_velo_proj->ksp, "velocity_gradient_projection_"));
112     {
113       PC pc;
114       PetscCall(KSPGetPC(grad_velo_proj->ksp, &pc));
115       PetscCall(PCSetType(pc, PCJACOBI));
116       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
117       PetscCall(KSPSetType(grad_velo_proj->ksp, KSPPREONLY));
118       // TODO Not sure if the option below are necessary
119       PetscCall(KSPSetConvergenceTest(grad_velo_proj->ksp, KSPConvergedSkip, NULL, NULL));
120     }
121     PetscCall(KSPSetOperators(grad_velo_proj->ksp, mat_mass, mat_mass));
122     PetscCall(KSPSetFromOptions(grad_velo_proj->ksp));
123   }
124 
125   PetscCallCeed(ceed, CeedBasisDestroy(&basis_grad_velo));
126   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
127   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_assemble));
128   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
129   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_assemble));
130   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
131   PetscFunctionReturn(PETSC_SUCCESS);
132 }
133 
134 PetscErrorCode VelocityGradientProjectionApply(User user, Vec Q_loc, Vec VelocityGradient) {
135   NodalProjectionData  grad_velo_proj = user->grad_velo_proj;
136   OperatorApplyContext l2_rhs_ctx     = grad_velo_proj->l2_rhs_ctx;
137 
138   PetscFunctionBeginUser;
139   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, VelocityGradient, l2_rhs_ctx));
140 
141   PetscCall(KSPSolve(grad_velo_proj->ksp, VelocityGradient, VelocityGradient));
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144