xref: /honee/src/smartsim/sgs_dd_training.c (revision 149fb5361f5198e41f87e8815a7e9dbfee84a96a)
1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2aa0b7f76SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3aa0b7f76SJames Wright //
4aa0b7f76SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5aa0b7f76SJames Wright //
6aa0b7f76SJames Wright // This file is part of CEED:  http://github.com/ceed
7aa0b7f76SJames Wright 
8aa0b7f76SJames Wright #include "../../qfunctions/sgs_dd_training.h"
9aa0b7f76SJames Wright 
10aa0b7f76SJames Wright #include <petscdmplex.h>
11aa0b7f76SJames Wright 
12*149fb536SJames Wright #include <navierstokes.h>
13*149fb536SJames Wright #include <smartsim.h>
14aa0b7f76SJames Wright 
15aa0b7f76SJames Wright typedef struct {
16aa0b7f76SJames Wright   CeedElemRestriction  elem_restr_grid_aniso;
17aa0b7f76SJames Wright   CeedVector           grid_aniso_ceed;
18aa0b7f76SJames Wright   CeedQFunctionContext sgs_dd_train_qfctx;
19aa0b7f76SJames Wright } *SGS_DD_TrainingSetupData;
20aa0b7f76SJames Wright 
21aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
22aa0b7f76SJames Wright   Ceed ceed;
23aa0b7f76SJames Wright 
24aa0b7f76SJames Wright   PetscFunctionBeginUser;
25aa0b7f76SJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_train_setup_data->elem_restr_grid_aniso, &ceed));
26aa0b7f76SJames Wright 
27aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_train_setup_data->elem_restr_grid_aniso));
28aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_train_setup_data->grid_aniso_ceed));
29aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_train_setup_data->sgs_dd_train_qfctx));
30aa0b7f76SJames Wright   PetscCall(PetscFree(sgs_dd_train_setup_data));
31aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
32aa0b7f76SJames Wright }
33aa0b7f76SJames Wright 
34aa0b7f76SJames Wright // @brief Create DM for storing data-drive SGS model inputs
35aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
36aa0b7f76SJames Wright   PetscSection section;
37aa0b7f76SJames Wright 
38aa0b7f76SJames Wright   PetscFunctionBeginUser;
39aa0b7f76SJames Wright   *num_components = 12;
40aa0b7f76SJames Wright 
41aa0b7f76SJames Wright   PetscCall(DMClone(dm_source, dm_dd_training));
42aa0b7f76SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data"));
43aa0b7f76SJames Wright 
44aa0b7f76SJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training));
45aa0b7f76SJames Wright 
46aa0b7f76SJames Wright   PetscCall(DMGetLocalSection(*dm_dd_training, &section));
47aa0b7f76SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data"));
48aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1"));
49aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2"));
50aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3"));
51aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4"));
52aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5"));
53aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6"));
54aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX"));
55aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY"));
56aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ"));
57aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ"));
58aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ"));
59aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY"));
60aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
61aa0b7f76SJames Wright };
62aa0b7f76SJames Wright 
63aa0b7f76SJames Wright // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes
64991aef52SJames Wright static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, User user, CeedData ceed_data, ProblemData problem,
65aa0b7f76SJames Wright                                                    SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
66aa0b7f76SJames Wright   SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train;
67aa0b7f76SJames Wright   CeedQFunction       qf_sgs_dd_train;
68aa0b7f76SJames Wright   CeedOperator        op_sgs_dd_train;
69aa0b7f76SJames Wright   CeedInt             num_comp_grad_velo, num_comp_grid_aniso;
70aa0b7f76SJames Wright   CeedVector          inv_multiplicity, filtered_fields;
71aa0b7f76SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train;
72aa0b7f76SJames Wright   DMLabel             domain_label = NULL;
73aa0b7f76SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
74aa0b7f76SJames Wright 
75aa0b7f76SJames Wright   PetscFunctionBeginUser;
76aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
77aa0b7f76SJames Wright 
78aa0b7f76SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, &elem_restr_sgs_train));
795930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, PETSC_TRUE,
805930f037SJames Wright                                    &elem_restr_inv_multiplicity, &inv_multiplicity));
81aa0b7f76SJames Wright 
82aa0b7f76SJames Wright   CeedElemRestriction elem_restr_filtered_state;
83aa0b7f76SJames Wright   CeedInt             num_comp_filtered_state;
84aa0b7f76SJames Wright   {  // -- Setup filtered velocity gradient projection
85aa0b7f76SJames Wright     CeedBasis         basis_filtered_state;
86aa0b7f76SJames Wright     CeedOperatorField op_field;
87aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->diff_filter->op_rhs_ctx->op, "v0", &op_field));
88aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_state));
89aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state));
90aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(op_field, &basis_filtered_state));
91aa0b7f76SJames Wright     PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state,
92aa0b7f76SJames Wright                                               &sgs_dd_train->filtered_grad_velo_proj));
93aa0b7f76SJames Wright     // Get velocity gradient information
94aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
95aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
96aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
97aa0b7f76SJames Wright   }
98aa0b7f76SJames Wright 
99aa0b7f76SJames Wright   CeedElemRestriction elem_restr_filtered_velo_prod;
100aa0b7f76SJames Wright   CeedInt             num_comp_filtered_velo_prod;
101aa0b7f76SJames Wright   {  // Get filtered velocity product information
102aa0b7f76SJames Wright     CeedOperatorField op_field;
103aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->diff_filter->op_rhs_ctx->op, "v1", &op_field));
104aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod));
105aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod));
106aa0b7f76SJames Wright   }
107aa0b7f76SJames Wright 
108aa0b7f76SJames Wright   // -- Create operator for generating training data at nodes
109aa0b7f76SJames Wright   // Differential Filter only provides filtered primitive variables
110aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim,
111aa0b7f76SJames Wright                                                   ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train));
112aa0b7f76SJames Wright 
113aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx));
114aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE));
115aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE));
116aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
117aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
118aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE));
119aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE));
120aa0b7f76SJames Wright 
121aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL));
122aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train));
123c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields));
124c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields));
125c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
126c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
127c6271fa9SJeremy L Thompson                                            sgs_dd_train_setup_data->grid_aniso_ceed));
128c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
129c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
130aa0b7f76SJames Wright 
131aa0b7f76SJames Wright   PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL,
132aa0b7f76SJames Wright                                        NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx));
133aa0b7f76SJames Wright 
134aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
135aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
136aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
137aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_train));
138aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_train));
139aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
140aa0b7f76SJames Wright }
141aa0b7f76SJames Wright 
142991aef52SJames Wright PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
143aa0b7f76SJames Wright   SGS_DDTrainingContext    sgsdd_train_qfctx;
144aa0b7f76SJames Wright   SGS_DD_TrainingSetupData sgs_dd_train_setup_data;
145aa0b7f76SJames Wright 
146aa0b7f76SJames Wright   PetscFunctionBeginUser;
147aa0b7f76SJames Wright   if (!user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
148aa0b7f76SJames Wright   if (!user->smartsim) PetscCall(SmartSimSetup(user));
149aa0b7f76SJames Wright 
150aa0b7f76SJames Wright   PetscCall(PetscNew(&sgsdd_train_qfctx));
151aa0b7f76SJames Wright   PetscCall(PetscNew(&sgs_dd_train_setup_data));
152aa0b7f76SJames Wright   PetscCall(PetscNew(&user->sgs_dd_train));
153aa0b7f76SJames Wright   SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train;
154aa0b7f76SJames Wright 
155aa0b7f76SJames Wright   sgs_dd_train->overwrite_training_data = PETSC_TRUE;
156aa0b7f76SJames Wright   sgs_dd_train->write_data_interval     = 1;
1576ea7c1aeSJames Wright   sgs_dd_train->num_filter_widths       = sizeof(sgs_dd_train->filter_widths) / sizeof(sgs_dd_train->filter_widths[0]);
158aa0b7f76SJames Wright   PetscOptionsBegin(user->comm, NULL, "SGS Data-Driven Training Options", NULL);
159aa0b7f76SJames Wright   PetscCall(PetscOptionsInt("-sgs_train_write_data_interval", "Number of timesteps between writing data into database", NULL,
160aa0b7f76SJames Wright                             sgs_dd_train->write_data_interval, &sgs_dd_train->write_data_interval, NULL));
161aa0b7f76SJames Wright   PetscCall(PetscOptionsBool("-sgs_train_overwrite_data", "Overwrite old training data in the database", NULL, sgs_dd_train->overwrite_training_data,
162aa0b7f76SJames Wright                              &sgs_dd_train->overwrite_training_data, NULL));
1636ea7c1aeSJames Wright   PetscCall(PetscOptionsRealArray("-sgs_train_filter_width_scales", "Scales of each filter width put into training database", NULL,
1646ea7c1aeSJames Wright                                   sgs_dd_train->filter_widths, &sgs_dd_train->num_filter_widths, NULL));
165aa0b7f76SJames Wright   PetscOptionsEnd();
166aa0b7f76SJames Wright 
167aa0b7f76SJames Wright   // -- Create DM for storing training data
168aa0b7f76SJames Wright   PetscCall(SGS_DD_TrainingCreateDM(user->dm, &sgs_dd_train->dm_dd_training, user->app_ctx->degree, user->app_ctx->q_extra,
169aa0b7f76SJames Wright                                     &sgs_dd_train->num_comp_dd_inputs));
170aa0b7f76SJames Wright 
171aa0b7f76SJames Wright   {  // -- Create QFunction Context
172aa0b7f76SJames Wright     NewtonianIdealGasContext gas;
173aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
174aa0b7f76SJames Wright     sgsdd_train_qfctx->gas = *gas;
175aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
176aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_train_setup_data->sgs_dd_train_qfctx));
177aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextSetData(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, CEED_USE_POINTER,
178aa0b7f76SJames Wright                                                     sizeof(*sgsdd_train_qfctx), sgsdd_train_qfctx));
179aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, FreeContextPetsc));
180aa0b7f76SJames Wright   }
181aa0b7f76SJames Wright 
182aa0b7f76SJames Wright   {  // -- Send training data array info to SmartRedis database
183aa0b7f76SJames Wright     PetscMPIInt  rank, num_ranks;
184aa0b7f76SJames Wright     SmartSimData smartsim = user->smartsim;
185aa0b7f76SJames Wright     PetscCallMPI(MPI_Comm_rank(user->comm, &rank));
186aa0b7f76SJames Wright     PetscCallMPI(MPI_Comm_size(user->comm, &num_ranks));
187aa0b7f76SJames Wright 
188aa0b7f76SJames Wright     {
189aa0b7f76SJames Wright       PetscSection global_section;
1907ff16c02SJames Wright       PetscInt     num_dofs, num_comps, local_min_max[2] = {0.}, global_min_max[2] = {0.};
1917ff16c02SJames Wright 
192aa0b7f76SJames Wright       PetscCall(DMGetGlobalSection(sgs_dd_train->dm_dd_training, &global_section));
193aa0b7f76SJames Wright       PetscCall(DMGetGlobalVectorInfo(sgs_dd_train->dm_dd_training, &num_dofs, NULL, NULL));
194aa0b7f76SJames Wright       PetscCall(PetscSectionGetFieldComponents(global_section, 0, &num_comps));
1957ff16c02SJames Wright       local_min_max[0] = num_dofs;
1967ff16c02SJames Wright       PetscCall(PetscGlobalMinMaxInt(user->comm, local_min_max, global_min_max));
1977ff16c02SJames Wright 
1987ff16c02SJames Wright       sgs_dd_train->training_data_array_dims[0] = global_min_max[0] / num_comps;
199aa0b7f76SJames Wright       sgs_dd_train->training_data_array_dims[1] = num_comps;
200aa0b7f76SJames Wright     }
201aa0b7f76SJames Wright 
202aa0b7f76SJames Wright     if (rank % smartsim->collocated_database_num_ranks == 0) {
2034fa1625aSJames Wright       {  // Communicate info on simulation size
2044fa1625aSJames Wright         const char tensor_name[]  = "sizeInfo";
205aa0b7f76SJames Wright         size_t     array_info_dim = 6;
206f4632befSRiccardo Balin         PetscInt64 array_info[6] = {0}, num_features = 6;
207aa0b7f76SJames Wright 
208aa0b7f76SJames Wright         array_info[0] = sgs_dd_train->training_data_array_dims[0];
209aa0b7f76SJames Wright         array_info[1] = sgs_dd_train->training_data_array_dims[1];
210aa0b7f76SJames Wright         array_info[2] = num_features;
211aa0b7f76SJames Wright         array_info[3] = num_ranks;
212aa0b7f76SJames Wright         array_info[4] = smartsim->collocated_database_num_ranks;
213aa0b7f76SJames Wright         array_info[5] = rank;
214aa0b7f76SJames Wright 
215ad2e713eSRiccardo Balin         PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
21643e9749fSJames Wright         PetscCallSmartRedis(
2174fa1625aSJames Wright             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), array_info, &array_info_dim, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
2184fa1625aSJames Wright         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
219ad2e713eSRiccardo Balin         PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
2204fa1625aSJames Wright       }
221aa0b7f76SJames Wright 
2224fa1625aSJames Wright       {  // Send array that communicates if tensors are overwritten in database
2234fa1625aSJames Wright         const char tensor_name[]       = "tensor-ow";
224f4632befSRiccardo Balin         PetscInt64 tensor_overwrite[2] = {sgs_dd_train->overwrite_training_data};
225aa0b7f76SJames Wright         size_t     dim_2[1]            = {2};
2264fa1625aSJames Wright 
227ad2e713eSRiccardo Balin         PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
22843e9749fSJames Wright         PetscCallSmartRedis(
2294fa1625aSJames Wright             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), tensor_overwrite, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
2304fa1625aSJames Wright         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
231ad2e713eSRiccardo Balin         PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
232aa0b7f76SJames Wright       }
2336ea7c1aeSJames Wright 
2346ea7c1aeSJames Wright       {  // Communicate number of filter widths used
2356ea7c1aeSJames Wright         const char tensor_name[]     = "num_filter_widths";
2366ea7c1aeSJames Wright         PetscInt64 num_filter_widths = sgs_dd_train->num_filter_widths;
2376ea7c1aeSJames Wright         size_t     dim_2             = 1;
2386ea7c1aeSJames Wright 
2396ea7c1aeSJames Wright         PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
24043e9749fSJames Wright         PetscCallSmartRedis(
2416ea7c1aeSJames Wright             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), &num_filter_widths, &dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
2426ea7c1aeSJames Wright         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
2436ea7c1aeSJames Wright         PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
2446ea7c1aeSJames Wright       }
245aa0b7f76SJames Wright     }
2464fa1625aSJames Wright   }
247aa0b7f76SJames Wright 
248aa0b7f76SJames Wright   // -- Compute and store anisotropy tensor
249aa0b7f76SJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_train_setup_data->elem_restr_grid_aniso,
250aa0b7f76SJames Wright                                                      &sgs_dd_train_setup_data->grid_aniso_ceed));
251aa0b7f76SJames Wright 
252aa0b7f76SJames Wright   // -- Create Nodal Evaluation Operator
253aa0b7f76SJames Wright   PetscCall(SetupTrainingDataCalculation(ceed, user, ceed_data, problem, sgs_dd_train_setup_data));
254aa0b7f76SJames Wright 
255aa0b7f76SJames Wright   PetscCall(SGS_DD_TrainingSetupDataDestroy(sgs_dd_train_setup_data));
256aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
257aa0b7f76SJames Wright }
258aa0b7f76SJames Wright 
259aa0b7f76SJames Wright PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx) {
260aa0b7f76SJames Wright   User                user         = (User)ctx;
261aa0b7f76SJames Wright   Ceed                ceed         = user->ceed;
262aa0b7f76SJames Wright   SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train;
263aa0b7f76SJames Wright   SmartSimData        smartsim     = user->smartsim;
264aa0b7f76SJames Wright   Vec                 TrainingData;
265ee81c423SRiccardo Balin   PetscMPIInt         rank;
266aa0b7f76SJames Wright 
267aa0b7f76SJames Wright   PetscFunctionBeginUser;
268ee81c423SRiccardo Balin 
269ee81c423SRiccardo Balin   PetscCallMPI(MPI_Comm_rank(user->comm, &rank));
270ee81c423SRiccardo Balin 
271aa0b7f76SJames Wright   if (step_num % sgs_dd_train->write_data_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
272aa0b7f76SJames Wright   PetscCall(DMGetGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData));
273aa0b7f76SJames Wright 
2746ea7c1aeSJames Wright   for (PetscInt filter_index = 0; filter_index < sgs_dd_train->num_filter_widths; filter_index++) {
275ad2e713eSRiccardo Balin     PetscCall(PetscLogEventBegin(FLUIDS_TrainDataCompute, 0, 0, 0, 0));
276aa0b7f76SJames Wright     {  // -- Compute and assemble training data
277aa0b7f76SJames Wright       Vec          FilteredVelocityGradient, FilteredFields, FilteredFields_loc;
278aa0b7f76SJames Wright       PetscMemType filtered_fields_mem_type;
279aa0b7f76SJames Wright       CeedVector   filtered_fields;
280aa0b7f76SJames Wright 
2816ea7c1aeSJames Wright       {  // Set filter width for the current solve
2826ea7c1aeSJames Wright         double       filter_width_scaling[3];
2836ea7c1aeSJames Wright         CeedOperator op_mat;
2846ea7c1aeSJames Wright         Mat          A_mat;
2856ea7c1aeSJames Wright 
2866ea7c1aeSJames Wright         for (int j = 0; j < 3; j++) filter_width_scaling[j] = sgs_dd_train->filter_widths[filter_index];
2876ea7c1aeSJames Wright         PetscCall(KSPGetOperators(user->diff_filter->ksp, &A_mat, NULL));
2886ea7c1aeSJames Wright         PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL));
2896ea7c1aeSJames Wright         PetscCall(CeedOperatorSetContextDouble(op_mat, user->diff_filter->filter_width_scaling_label, filter_width_scaling));
2906ea7c1aeSJames Wright       }
2916ea7c1aeSJames Wright 
292aa0b7f76SJames Wright       PetscCall(DMGetGlobalVector(user->diff_filter->dm_filter, &FilteredFields));
293aa0b7f76SJames Wright       PetscCall(DMGetLocalVector(user->diff_filter->dm_filter, &FilteredFields_loc));
294aa0b7f76SJames Wright 
295aa0b7f76SJames Wright       PetscCall(DifferentialFilterApply(user, solution_time, Q, FilteredFields));
296aa0b7f76SJames Wright       PetscCall(DMGlobalToLocal(user->diff_filter->dm_filter, FilteredFields, INSERT_VALUES, FilteredFields_loc));
297aa0b7f76SJames Wright 
298aa0b7f76SJames Wright       PetscCall(DMGetGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient));
299aa0b7f76SJames Wright       PetscCall(VelocityGradientProjectionApply(sgs_dd_train->filtered_grad_velo_proj, FilteredFields_loc, FilteredVelocityGradient));
300aa0b7f76SJames Wright 
301aa0b7f76SJames Wright       {
302aa0b7f76SJames Wright         CeedOperatorField op_field;
3036ea7c1aeSJames Wright 
304aa0b7f76SJames Wright         PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->op_training_data_calc_ctx->op, "q", &op_field));
305aa0b7f76SJames Wright         PetscCallCeed(ceed, CeedOperatorFieldGetVector(op_field, &filtered_fields));
306aa0b7f76SJames Wright       }
3076ea7c1aeSJames Wright 
308a7dac1d5SJames Wright       PetscCall(VecPetscToCeed(FilteredFields_loc, &filtered_fields_mem_type, filtered_fields));  // filtered_fields is an implicit input
309aa0b7f76SJames Wright       PetscCall(ApplyCeedOperatorGlobalToGlobal(FilteredVelocityGradient, TrainingData, sgs_dd_train->op_training_data_calc_ctx));
310a7dac1d5SJames Wright       PetscCall(VecCeedToPetsc(filtered_fields, filtered_fields_mem_type, FilteredFields_loc));
311aa0b7f76SJames Wright 
312aa0b7f76SJames Wright       PetscCall(DMRestoreGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient));
313aa0b7f76SJames Wright       PetscCall(DMRestoreGlobalVector(user->diff_filter->dm_filter, &FilteredFields));
314aa0b7f76SJames Wright       PetscCall(DMRestoreLocalVector(user->diff_filter->dm_filter, &FilteredFields_loc));
315aa0b7f76SJames Wright     }
316ad2e713eSRiccardo Balin     PetscCall(PetscLogEventEnd(FLUIDS_TrainDataCompute, 0, 0, 0, 0));
317aa0b7f76SJames Wright 
318aa0b7f76SJames Wright     {  // -- Send training data to SmartSim
319aa0b7f76SJames Wright       char   array_key[PETSC_MAX_PATH_LEN];
320aa0b7f76SJames Wright       size_t array_key_len;
321aa0b7f76SJames Wright 
322aa0b7f76SJames Wright       if (sgs_dd_train->overwrite_training_data) {
3236ea7c1aeSJames Wright         PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index));
324aa0b7f76SJames Wright       } else {
3256ea7c1aeSJames Wright         PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index));
326aa0b7f76SJames Wright       }
327aa0b7f76SJames Wright       PetscCall(PetscStrlen(array_key, &array_key_len));
328aa0b7f76SJames Wright 
329aa0b7f76SJames Wright       {
330aa0b7f76SJames Wright         const PetscScalar *training_data;
331aa0b7f76SJames Wright         PetscCall(VecGetArrayRead(TrainingData, &training_data));
332ad2e713eSRiccardo Balin         PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Train, 0, 0, 0, 0));
33343e9749fSJames Wright         PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2,
334aa0b7f76SJames Wright                                        SRTensorTypeDouble, SRMemLayoutContiguous));
335ad2e713eSRiccardo Balin         PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Train, 0, 0, 0, 0));
336aa0b7f76SJames Wright         PetscCall(VecRestoreArrayRead(TrainingData, &training_data));
337aa0b7f76SJames Wright       }
338ee81c423SRiccardo Balin     }
339ee81c423SRiccardo Balin   }
340aa0b7f76SJames Wright 
341aa0b7f76SJames Wright   if (rank % smartsim->collocated_database_num_ranks == 0) {
3426ea7c1aeSJames Wright     const char tensor_name[] = "step";
343aa0b7f76SJames Wright     size_t     dim_2[1]      = {2};
344f4632befSRiccardo Balin     PetscInt64 step_array[2] = {step_num, step_num};
3456ea7c1aeSJames Wright 
346ad2e713eSRiccardo Balin     PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
34743e9749fSJames Wright     PetscCallSmartRedis(
3486ea7c1aeSJames Wright         put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
349ad2e713eSRiccardo Balin     PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
350aa0b7f76SJames Wright   }
351aa0b7f76SJames Wright 
352aa0b7f76SJames Wright   PetscCall(DMRestoreGlobalVector(user->sgs_dd_train->dm_dd_training, &TrainingData));
353aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
354aa0b7f76SJames Wright }
355aa0b7f76SJames Wright 
356632a41e1SJames Wright PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) {
357632a41e1SJames Wright   User         user;
358632a41e1SJames Wright   const char   check_run_key[]   = "check-run";
359632a41e1SJames Wright   PetscReal    check_run[2]      = {1};
360632a41e1SJames Wright   const size_t check_run_dims[1] = {2};
361632a41e1SJames Wright   size_t       check_run_key_size;
362632a41e1SJames Wright 
363632a41e1SJames Wright   PetscFunctionBeginUser;
364632a41e1SJames Wright   PetscCall(PetscStrlen(check_run_key, &check_run_key_size));
365632a41e1SJames Wright   PetscCall(TSGetApplicationContext(ts, &user));
366632a41e1SJames Wright   SmartSimData smartsim = user->smartsim;
367632a41e1SJames Wright 
368ad2e713eSRiccardo Balin   PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
36943e9749fSJames Wright   PetscCallSmartRedis(
370632a41e1SJames Wright       unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous));
371ad2e713eSRiccardo Balin   PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0));
372632a41e1SJames Wright   if (check_run[0] == 0) {
373632a41e1SJames Wright     PetscCall(PetscPrintf(user->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n"));
374632a41e1SJames Wright     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
375632a41e1SJames Wright   }
376632a41e1SJames Wright 
377632a41e1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
378632a41e1SJames Wright }
379632a41e1SJames Wright 
380aa0b7f76SJames Wright PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train) {
381aa0b7f76SJames Wright   PetscFunctionBeginUser;
382aa0b7f76SJames Wright   if (!sgs_dd_train) PetscFunctionReturn(PETSC_SUCCESS);
383aa0b7f76SJames Wright 
384aa0b7f76SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_train->op_training_data_calc_ctx));
385aa0b7f76SJames Wright   PetscCall(NodalProjectionDataDestroy(sgs_dd_train->filtered_grad_velo_proj));
386aa0b7f76SJames Wright   PetscCall(DMDestroy(&sgs_dd_train->dm_dd_training));
387aa0b7f76SJames Wright   PetscCall(PetscFree(sgs_dd_train));
388aa0b7f76SJames Wright 
389aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
390aa0b7f76SJames Wright }
391