xref: /honee/problems/sgs_dd_model.c (revision ee1455b73fdf3c948d552fc0a19cf21e2ac18abb)
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 
8 #include "../qfunctions/sgs_dd_model.h"
9 
10 #include <petscdmplex.h>
11 
12 #include "../navierstokes.h"
13 
14 typedef struct {
15   CeedElemRestriction  elem_restr_grid_aniso, elem_restr_sgs;
16   CeedVector           grid_aniso_ceed;
17   CeedQFunctionContext sgsdd_qfctx;
18 } *SGS_DD_ModelSetupData;
19 
20 PetscErrorCode SGS_DD_ModelSetupDataDestroy(SGS_DD_ModelSetupData sgs_dd_setup_data) {
21   PetscFunctionBeginUser;
22   CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso);
23   CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs);
24   CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed);
25   CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx);
26 
27   PetscCall(PetscFree(sgs_dd_setup_data));
28   PetscFunctionReturn(0);
29 }
30 
31 // @brief Create DM for storing subgrid stress at nodes
32 PetscErrorCode SGS_DD_ModelCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt *num_components) {
33   PetscFE      fe;
34   PetscSection section;
35   PetscInt     dim;
36 
37   PetscFunctionBeginUser;
38   *num_components = 6;
39 
40   PetscCall(DMClone(dm_source, dm_sgs));
41   PetscCall(DMGetDimension(*dm_sgs, &dim));
42   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
43 
44   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, *num_components, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
45   PetscCall(PetscObjectSetName((PetscObject)fe, "Subgrid Stress Projection"));
46   PetscCall(DMAddField(*dm_sgs, NULL, (PetscObject)fe));
47   PetscCall(DMCreateDS(*dm_sgs));
48   PetscCall(DMPlexSetClosurePermutationTensor(*dm_sgs, PETSC_DETERMINE, NULL));
49 
50   PetscCall(DMGetLocalSection(*dm_sgs, &section));
51   PetscCall(PetscSectionSetFieldName(section, 0, ""));
52   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
53   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
54   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
55   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
56   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
57   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
58 
59   PetscCall(PetscFEDestroy(&fe));
60 
61   PetscFunctionReturn(0);
62 };
63 
64 // @brief Create CeedOperator to calculate data-drive SGS at nodes
65 PetscErrorCode SGS_DD_ModelSetupNodalEvaluation(Ceed ceed, User user, CeedData ceed_data, SGS_DD_ModelSetupData sgs_dd_setup_data) {
66   SGS_DD_Data         sgs_dd_data = user->sgs_dd_data;
67   CeedQFunction       qf_multiplicity, qf_sgs_dd_nodal;
68   CeedOperator        op_multiplicity, op_sgs_dd_nodal;
69   CeedInt             num_elem, elem_size, num_comp_q, dim, num_qpts_1d, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
70   CeedVector          multiplicity, inv_multiplicity, grad_velo_ceed;
71   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
72 
73   PetscFunctionBeginUser;
74   PetscCall(DMGetDimension(user->dm, &dim));
75   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x);
76   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
77   CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso);
78   CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &num_elem);
79   CeedElemRestrictionGetElementSize(ceed_data->elem_restr_q, &elem_size);
80   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &num_qpts_1d);
81 
82   {  // Get velocity gradient information
83     CeedOperatorField op_field;
84     CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field);
85     CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo);
86     CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo);
87     CeedElemRestrictionCreateVector(elem_restr_grad_velo, &grad_velo_ceed, NULL);
88   }
89 
90   PetscCall(GetRestrictionForDomain(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, num_qpts_1d, 0, &elem_restr_sgs, NULL, NULL));
91   CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL);
92 
93   // -- Create inverse multiplicity for correcting nodal assembly
94   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &multiplicity, NULL);
95   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, multiplicity);
96   CeedElemRestrictionCreateStrided(ceed, num_elem, elem_size, 1, num_elem * elem_size, CEED_STRIDES_BACKEND, &elem_restr_inv_multiplicity);
97   CeedElemRestrictionCreateVector(elem_restr_inv_multiplicity, &inv_multiplicity, NULL);
98 
99   CeedQFunctionCreateInterior(ceed, 1, InverseMultiplicity, InverseMultiplicity_loc, &qf_multiplicity);
100   CeedQFunctionAddInput(qf_multiplicity, "multiplicity", num_comp_q, CEED_EVAL_NONE);
101   CeedQFunctionAddOutput(qf_multiplicity, "inverse multiplicity", 1, CEED_EVAL_NONE);
102 
103   CeedOperatorCreate(ceed, qf_multiplicity, NULL, NULL, &op_multiplicity);
104   CeedOperatorSetName(op_multiplicity, "SGS DD Model - Create Multiplicity Scaling");
105   CeedOperatorSetField(op_multiplicity, "multiplicity", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
106   CeedOperatorSetField(op_multiplicity, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
107   CeedOperatorSetNumQuadraturePoints(op_multiplicity, elem_size);
108 
109   CeedOperatorApply(op_multiplicity, multiplicity, inv_multiplicity, CEED_REQUEST_IMMEDIATE);
110 
111   // -- Create operator for SGS DD model nodal evaluation
112   switch (user->phys->state_var) {
113     case STATEVAR_PRIMITIVE:
114       CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicNodal_Prim, ComputeSGS_DDAnisotropicNodal_Prim_loc, &qf_sgs_dd_nodal);
115       break;
116     case STATEVAR_CONSERVATIVE:
117       CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicNodal_Conserv, ComputeSGS_DDAnisotropicNodal_Conserv_loc, &qf_sgs_dd_nodal);
118       break;
119     default:
120       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP,
121               "Anisotropic data-driven SGS nodal evaluation not available for chosen state variable");
122   }
123 
124   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
125   CeedBasis basis_x_to_q;
126   PetscCall(CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q));
127 
128   CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx);
129   CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE);
130   CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP);
131   CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE);
132   CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE);
133   CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE);
134   CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE);
135 
136   CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal);
137   CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, user->q_ceed);
138   CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord);
139   CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
140   CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_COLLOCATED,
141                        sgs_dd_setup_data->grid_aniso_ceed);
142   CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_COLLOCATED, inv_multiplicity);
143   CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
144 
145   PetscCall(OperatorApplyContextCreate(user->grad_velo_proj->dm, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_nodal, grad_velo_ceed,
146                                        sgs_dd_data->sgs_nodal_ceed, NULL, NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
147 
148   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
149 
150   CeedVectorDestroy(&multiplicity);
151   CeedVectorDestroy(&inv_multiplicity);
152   CeedVectorDestroy(&grad_velo_ceed);
153   CeedBasisDestroy(&basis_x_to_q);
154   CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity);
155   CeedQFunctionDestroy(&qf_multiplicity);
156   CeedQFunctionDestroy(&qf_sgs_dd_nodal);
157   CeedOperatorDestroy(&op_multiplicity);
158   CeedOperatorDestroy(&op_sgs_dd_nodal);
159   PetscFunctionReturn(0);
160 }
161 
162 // @brief B = A^T, A is NxM, B is MxN
163 PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
164   PetscFunctionBeginUser;
165   for (PetscInt i = 0; i < N; i++) {
166     for (PetscInt j = 0; j < M; j++) {
167       B[j * N + i] = A[i * M + j];
168     }
169   }
170   PetscFunctionReturn(0);
171 }
172 
173 // @brief Read neural network coefficients from file and put into context struct
174 PetscErrorCode SGS_DD_ModelContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SGS_DDModelContext *psgsdd_ctx) {
175   SGS_DDModelContext sgsdd_ctx;
176   PetscInt           num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
177   char               file_path[PETSC_MAX_PATH_LEN];
178   PetscScalar       *temp;
179 
180   PetscFunctionBeginUser;
181   {
182     SGS_DDModelContext sgsdd_temp;
183     PetscCall(PetscNew(&sgsdd_temp));
184     *sgsdd_temp                     = **psgsdd_ctx;
185     sgsdd_temp->offsets.bias1       = 0;
186     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
187     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
188     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
189     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
190     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
191     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
192     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
193     *sgsdd_ctx = *sgsdd_temp;
194     PetscCall(PetscFree(sgsdd_temp));
195   }
196 
197   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
198   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
199   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
200   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
201   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
202   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
203 
204   {
205     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
206     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
207     PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, temp));
208     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
209     PetscCall(PetscFree(temp));
210   }
211   {
212     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
213     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
214     PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, temp));
215     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
216     PetscCall(PetscFree(temp));
217   }
218 
219   PetscCall(PetscFree(*psgsdd_ctx));
220   *psgsdd_ctx = sgsdd_ctx;
221   PetscFunctionReturn(0);
222 }
223 
224 PetscErrorCode SGS_DD_ModelSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
225   PetscReal                alpha = 0;
226   SGS_DDModelContext       sgsdd_ctx;
227   MPI_Comm                 comm                           = user->comm;
228   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
229   SGS_DD_ModelSetupData    sgs_dd_setup_data;
230   NewtonianIdealGasContext gas;
231   PetscFunctionBeginUser;
232 
233   PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem));
234 
235   PetscCall(PetscNew(&sgsdd_ctx));
236 
237   PetscOptionsBegin(comm, NULL, "SGS Data-Drive Model Options", NULL);
238   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
239   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
240                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
241   PetscOptionsEnd();
242 
243   sgsdd_ctx->num_layers  = 2;
244   sgsdd_ctx->num_inputs  = 6;
245   sgsdd_ctx->num_outputs = 6;
246   sgsdd_ctx->num_neurons = 20;
247   sgsdd_ctx->alpha       = alpha;
248 
249   // PetscCall(SGS_DD_ModelContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
250 
251   // -- Create DM for storing SGS tensor at nodes
252   PetscCall(PetscNew(&user->sgs_dd_data));
253   PetscCall(SGS_DD_ModelCreateDM(user->dm, &user->sgs_dd_data->dm_sgs, user->app_ctx->degree, &user->sgs_dd_data->num_comp_sgs));
254 
255   PetscCall(PetscNew(&sgs_dd_setup_data));
256 
257   CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas);
258   sgsdd_ctx->gas = *gas;
259   CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas);
260   CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx);
261   CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx);
262   CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc);
263 
264   // -- Compute and store anisotropy tensor
265   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso,
266                                                      &sgs_dd_setup_data->grid_aniso_ceed));
267 
268   // -- Create Nodal Evaluation Operator
269   PetscCall(SGS_DD_ModelSetupNodalEvaluation(ceed, user, ceed_data, sgs_dd_setup_data));
270 
271   PetscCall(SGS_DD_ModelSetupDataDestroy(sgs_dd_setup_data));
272   PetscFunctionReturn(0);
273 }
274 
275 PetscErrorCode SGS_DD_DataDestroy(SGS_DD_Data sgs_dd_data) {
276   PetscFunctionBeginUser;
277   if (!sgs_dd_data) PetscFunctionReturn(0);
278 
279   CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed);
280   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
281   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
282   PetscCall(PetscFree(sgs_dd_data));
283 
284   PetscFunctionReturn(0);
285 }
286