xref: /honee/src/smartsim/sgs_dd_training.c (revision cb8a476cda62cb3d4c49be98bd55cb366864d7e3)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3aa0b7f76SJames Wright 
4aa0b7f76SJames Wright #include "../../qfunctions/sgs_dd_training.h"
5aa0b7f76SJames Wright 
6aa0b7f76SJames Wright #include <petscdmplex.h>
7aa0b7f76SJames Wright 
8149fb536SJames Wright #include <navierstokes.h>
9149fb536SJames Wright #include <smartsim.h>
10aa0b7f76SJames Wright 
11aa0b7f76SJames Wright typedef struct {
1239169b57SJames Wright   DM                   dm_dd_training;
1339169b57SJames Wright   PetscInt             num_comp_dd_inputs, write_data_interval, num_filter_widths;
1439169b57SJames Wright   PetscScalar          filter_widths[16];
1539169b57SJames Wright   OperatorApplyContext op_training_data_calc_ctx;
16*cb8a476cSJames Wright   DiffFilterData       diff_filter;
1739169b57SJames Wright   NodalProjectionData  filtered_grad_velo_proj;
1839169b57SJames Wright   size_t               training_data_array_dims[2];
1939169b57SJames Wright   PetscBool            overwrite_training_data;
2039169b57SJames Wright } *SGS_DD_TrainingData;
2139169b57SJames Wright 
2239169b57SJames Wright #define SGS_DD_TRAIN_KEY "SGS Data Driven Training"
2339169b57SJames Wright 
2439169b57SJames Wright PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData *sgs_dd_train) {
2539169b57SJames Wright   SGS_DD_TrainingData sgs_dd_train_ = *sgs_dd_train;
2639169b57SJames Wright 
2739169b57SJames Wright   PetscFunctionBeginUser;
2839169b57SJames Wright   if (!sgs_dd_train_) PetscFunctionReturn(PETSC_SUCCESS);
2939169b57SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_train_->op_training_data_calc_ctx));
3039169b57SJames Wright   PetscCall(NodalProjectionDataDestroy(&sgs_dd_train_->filtered_grad_velo_proj));
3139169b57SJames Wright   PetscCall(DMDestroy(&sgs_dd_train_->dm_dd_training));
32*cb8a476cSJames Wright   PetscCall(DifferentialFilterDataDestroy(&sgs_dd_train_->diff_filter));
3339169b57SJames Wright   PetscCall(PetscFree(sgs_dd_train));
3439169b57SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3539169b57SJames Wright }
3639169b57SJames Wright 
3739169b57SJames Wright typedef struct {
38aa0b7f76SJames Wright   CeedElemRestriction  elem_restr_grid_aniso;
39aa0b7f76SJames Wright   CeedVector           grid_aniso_ceed;
40aa0b7f76SJames Wright   CeedQFunctionContext sgs_dd_train_qfctx;
41aa0b7f76SJames Wright } *SGS_DD_TrainingSetupData;
42aa0b7f76SJames Wright 
43aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
44aa0b7f76SJames Wright   Ceed ceed;
45aa0b7f76SJames Wright 
46aa0b7f76SJames Wright   PetscFunctionBeginUser;
47aa0b7f76SJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_train_setup_data->elem_restr_grid_aniso, &ceed));
48aa0b7f76SJames Wright 
49aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_train_setup_data->elem_restr_grid_aniso));
50aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_train_setup_data->grid_aniso_ceed));
51aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_train_setup_data->sgs_dd_train_qfctx));
52aa0b7f76SJames Wright   PetscCall(PetscFree(sgs_dd_train_setup_data));
53519781aeSJames Wright   PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed");
54aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
55aa0b7f76SJames Wright }
56aa0b7f76SJames Wright 
57aa0b7f76SJames Wright // @brief Create DM for storing data-drive SGS model inputs
58aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
59aa0b7f76SJames Wright   PetscSection section;
60aa0b7f76SJames Wright 
61aa0b7f76SJames Wright   PetscFunctionBeginUser;
62aa0b7f76SJames Wright   *num_components = 12;
63aa0b7f76SJames Wright 
64aa0b7f76SJames Wright   PetscCall(DMClone(dm_source, dm_dd_training));
650dee9b8eSJames Wright   PetscCall(DMSetMatrixPreallocateSkip(*dm_dd_training, PETSC_TRUE));
66aa0b7f76SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data"));
67aa0b7f76SJames Wright 
68aa0b7f76SJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training));
69aa0b7f76SJames Wright 
70aa0b7f76SJames Wright   PetscCall(DMGetLocalSection(*dm_dd_training, &section));
71aa0b7f76SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data"));
72aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1"));
73aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2"));
74aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3"));
75aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4"));
76aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5"));
77aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6"));
78aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX"));
79aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY"));
80aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ"));
81aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ"));
82aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ"));
83aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY"));
84aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
85aa0b7f76SJames Wright };
86aa0b7f76SJames Wright 
87aa0b7f76SJames Wright // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes
88e3663b90SJames Wright static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, Honee honee, ProblemData problem, SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
8939169b57SJames Wright   SGS_DD_TrainingData sgs_dd_train;
90aa0b7f76SJames Wright   CeedQFunction       qf_sgs_dd_train;
91aa0b7f76SJames Wright   CeedOperator        op_sgs_dd_train;
92aa0b7f76SJames Wright   CeedInt             num_comp_grad_velo, num_comp_grid_aniso;
93aa0b7f76SJames Wright   CeedVector          inv_multiplicity, filtered_fields;
94aa0b7f76SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train;
95aa0b7f76SJames Wright   DMLabel             domain_label = NULL;
96aa0b7f76SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
97aa0b7f76SJames Wright 
98aa0b7f76SJames Wright   PetscFunctionBeginUser;
9939169b57SJames Wright   PetscCall(PetscObjectContainerQuery((PetscObject)honee, SGS_DD_TRAIN_KEY, &sgs_dd_train));
100aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
101aa0b7f76SJames Wright 
102aa0b7f76SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, &elem_restr_sgs_train));
1035930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, PETSC_TRUE,
1045930f037SJames Wright                                    &elem_restr_inv_multiplicity, &inv_multiplicity));
105aa0b7f76SJames Wright 
106aa0b7f76SJames Wright   CeedElemRestriction elem_restr_filtered_state;
107aa0b7f76SJames Wright   CeedInt             num_comp_filtered_state;
108aa0b7f76SJames Wright   {  // -- Setup filtered velocity gradient projection
109aa0b7f76SJames Wright     CeedBasis         basis_filtered_state;
110aa0b7f76SJames Wright     CeedOperatorField op_field;
111*cb8a476cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->diff_filter->op_rhs_ctx->op, "v0", &op_field));
11201e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_filtered_state, &basis_filtered_state, NULL));
113aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state));
114e3663b90SJames Wright     PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state,
115aa0b7f76SJames Wright                                               &sgs_dd_train->filtered_grad_velo_proj));
116fff85bd3SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_filtered_state));
117aa0b7f76SJames Wright     // Get velocity gradient information
118aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
119aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
120aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
121aa0b7f76SJames Wright   }
122aa0b7f76SJames Wright 
123aa0b7f76SJames Wright   CeedElemRestriction elem_restr_filtered_velo_prod;
124aa0b7f76SJames Wright   CeedInt             num_comp_filtered_velo_prod;
125aa0b7f76SJames Wright   {  // Get filtered velocity product information
126aa0b7f76SJames Wright     CeedOperatorField op_field;
127*cb8a476cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->diff_filter->op_rhs_ctx->op, "v1", &op_field));
128aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod));
129aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod));
130aa0b7f76SJames Wright   }
131aa0b7f76SJames Wright 
132aa0b7f76SJames Wright   // -- Create operator for generating training data at nodes
133aa0b7f76SJames Wright   // Differential Filter only provides filtered primitive variables
134aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim,
135aa0b7f76SJames Wright                                                   ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train));
136aa0b7f76SJames Wright 
137aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx));
138aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE));
139aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE));
140aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
141aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
142aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE));
143aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE));
144aa0b7f76SJames Wright 
145aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL));
146aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train));
147c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields));
148c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields));
149c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
150c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
151c6271fa9SJeremy L Thompson                                            sgs_dd_train_setup_data->grid_aniso_ceed));
152c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
153c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
154aa0b7f76SJames Wright 
155aa0b7f76SJames Wright   PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL,
156aa0b7f76SJames Wright                                        NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx));
157aa0b7f76SJames Wright 
158aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
159aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
160aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
161fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_state));
162fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
163fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_velo_prod));
164aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_train));
165aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_train));
166aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
167aa0b7f76SJames Wright }
168aa0b7f76SJames Wright 
169e3663b90SJames Wright PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem) {
1707f3a2123SJames Wright   SGS_DDTrainingContext    sgsdd_train_ctx;
171aa0b7f76SJames Wright   SGS_DD_TrainingSetupData sgs_dd_train_setup_data;
17239169b57SJames Wright   SGS_DD_TrainingData      sgs_dd_train;
173aa0b7f76SJames Wright 
174aa0b7f76SJames Wright   PetscFunctionBeginUser;
1750c373b74SJames Wright   if (!honee->smartsim) PetscCall(SmartSimSetup(honee));
176aa0b7f76SJames Wright 
1777f3a2123SJames Wright   PetscCall(PetscNew(&sgsdd_train_ctx));
178aa0b7f76SJames Wright   PetscCall(PetscNew(&sgs_dd_train_setup_data));
17939169b57SJames Wright   PetscCall(PetscNew(&sgs_dd_train));
18039169b57SJames Wright   PetscCall(PetscObjectContainerCompose((PetscObject)honee, SGS_DD_TRAIN_KEY, sgs_dd_train, (PetscCtxDestroyFn *)SGS_DD_TrainingDataDestroy));
181*cb8a476cSJames Wright   PetscCall(DifferentialFilterSetup(honee, &sgs_dd_train->diff_filter));
182aa0b7f76SJames Wright 
183aa0b7f76SJames Wright   sgs_dd_train->overwrite_training_data = PETSC_TRUE;
184aa0b7f76SJames Wright   sgs_dd_train->write_data_interval     = 1;
1856ea7c1aeSJames Wright   sgs_dd_train->num_filter_widths       = sizeof(sgs_dd_train->filter_widths) / sizeof(sgs_dd_train->filter_widths[0]);
1860c373b74SJames Wright   PetscOptionsBegin(honee->comm, NULL, "SGS Data-Driven Training Options", NULL);
187aa0b7f76SJames Wright   PetscCall(PetscOptionsInt("-sgs_train_write_data_interval", "Number of timesteps between writing data into database", NULL,
188aa0b7f76SJames Wright                             sgs_dd_train->write_data_interval, &sgs_dd_train->write_data_interval, NULL));
189aa0b7f76SJames Wright   PetscCall(PetscOptionsBool("-sgs_train_overwrite_data", "Overwrite old training data in the database", NULL, sgs_dd_train->overwrite_training_data,
190aa0b7f76SJames Wright                              &sgs_dd_train->overwrite_training_data, NULL));
1916ea7c1aeSJames Wright   PetscCall(PetscOptionsRealArray("-sgs_train_filter_width_scales", "Scales of each filter width put into training database", NULL,
1926ea7c1aeSJames Wright                                   sgs_dd_train->filter_widths, &sgs_dd_train->num_filter_widths, NULL));
193aa0b7f76SJames Wright   PetscOptionsEnd();
194aa0b7f76SJames Wright 
195aa0b7f76SJames Wright   // -- Create DM for storing training data
1960c373b74SJames Wright   PetscCall(SGS_DD_TrainingCreateDM(honee->dm, &sgs_dd_train->dm_dd_training, honee->app_ctx->degree, honee->app_ctx->q_extra,
197aa0b7f76SJames Wright                                     &sgs_dd_train->num_comp_dd_inputs));
198aa0b7f76SJames Wright 
199aa0b7f76SJames Wright   {  // -- Create QFunction Context
200aa0b7f76SJames Wright     NewtonianIdealGasContext gas;
201e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas));
2027f3a2123SJames Wright     sgsdd_train_ctx->gas = *gas;
203e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas));
2040c373b74SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_train_setup_data->sgs_dd_train_qfctx));
205aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextSetData(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, CEED_USE_POINTER,
2067f3a2123SJames Wright                                                     sizeof(*sgsdd_train_ctx), sgsdd_train_ctx));
207aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, FreeContextPetsc));
208aa0b7f76SJames Wright   }
209aa0b7f76SJames Wright 
210aa0b7f76SJames Wright   {  // -- Send training data array info to SmartRedis database
211aa0b7f76SJames Wright     PetscMPIInt  rank, num_ranks;
2120c373b74SJames Wright     SmartSimData smartsim = honee->smartsim;
2130c373b74SJames Wright     PetscCallMPI(MPI_Comm_rank(honee->comm, &rank));
2140c373b74SJames Wright     PetscCallMPI(MPI_Comm_size(honee->comm, &num_ranks));
215aa0b7f76SJames Wright 
216aa0b7f76SJames Wright     {
217aa0b7f76SJames Wright       PetscSection global_section;
2187ff16c02SJames Wright       PetscInt     num_dofs, num_comps, local_min_max[2] = {0.}, global_min_max[2] = {0.};
2197ff16c02SJames Wright 
220aa0b7f76SJames Wright       PetscCall(DMGetGlobalSection(sgs_dd_train->dm_dd_training, &global_section));
221aa0b7f76SJames Wright       PetscCall(DMGetGlobalVectorInfo(sgs_dd_train->dm_dd_training, &num_dofs, NULL, NULL));
222aa0b7f76SJames Wright       PetscCall(PetscSectionGetFieldComponents(global_section, 0, &num_comps));
2237ff16c02SJames Wright       local_min_max[0] = num_dofs;
2240c373b74SJames Wright       PetscCall(PetscGlobalMinMaxInt(honee->comm, local_min_max, global_min_max));
2257ff16c02SJames Wright 
2267ff16c02SJames Wright       sgs_dd_train->training_data_array_dims[0] = global_min_max[0] / num_comps;
227aa0b7f76SJames Wright       sgs_dd_train->training_data_array_dims[1] = num_comps;
228aa0b7f76SJames Wright     }
229aa0b7f76SJames Wright 
230aa0b7f76SJames Wright     if (rank % smartsim->collocated_database_num_ranks == 0) {
2314fa1625aSJames Wright       {  // Communicate info on simulation size
2324fa1625aSJames Wright         const char tensor_name[]  = "sizeInfo";
233aa0b7f76SJames Wright         size_t     array_info_dim = 6;
234f4632befSRiccardo Balin         PetscInt64 array_info[6] = {0}, num_features = 6;
235aa0b7f76SJames Wright 
236aa0b7f76SJames Wright         array_info[0] = sgs_dd_train->training_data_array_dims[0];
237aa0b7f76SJames Wright         array_info[1] = sgs_dd_train->training_data_array_dims[1];
238aa0b7f76SJames Wright         array_info[2] = num_features;
239aa0b7f76SJames Wright         array_info[3] = num_ranks;
240aa0b7f76SJames Wright         array_info[4] = smartsim->collocated_database_num_ranks;
241aa0b7f76SJames Wright         array_info[5] = rank;
242aa0b7f76SJames Wright 
243ea615d4cSJames Wright         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
24443e9749fSJames Wright         PetscCallSmartRedis(
2454fa1625aSJames Wright             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), array_info, &array_info_dim, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
2464fa1625aSJames Wright         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
247ea615d4cSJames Wright         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
2484fa1625aSJames Wright       }
249aa0b7f76SJames Wright 
2504fa1625aSJames Wright       {  // Send array that communicates if tensors are overwritten in database
2514fa1625aSJames Wright         const char tensor_name[]       = "tensor-ow";
252f4632befSRiccardo Balin         PetscInt64 tensor_overwrite[2] = {sgs_dd_train->overwrite_training_data};
253aa0b7f76SJames Wright         size_t     dim_2[1]            = {2};
2544fa1625aSJames Wright 
255ea615d4cSJames Wright         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
25643e9749fSJames Wright         PetscCallSmartRedis(
2574fa1625aSJames Wright             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), tensor_overwrite, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
2584fa1625aSJames Wright         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
259ea615d4cSJames Wright         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
260aa0b7f76SJames Wright       }
2616ea7c1aeSJames Wright 
2626ea7c1aeSJames Wright       {  // Communicate number of filter widths used
2636ea7c1aeSJames Wright         const char tensor_name[]     = "num_filter_widths";
2646ea7c1aeSJames Wright         PetscInt64 num_filter_widths = sgs_dd_train->num_filter_widths;
2656ea7c1aeSJames Wright         size_t     dim_2             = 1;
2666ea7c1aeSJames Wright 
267ea615d4cSJames Wright         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
26843e9749fSJames Wright         PetscCallSmartRedis(
2696ea7c1aeSJames Wright             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), &num_filter_widths, &dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
2706ea7c1aeSJames Wright         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
271ea615d4cSJames Wright         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
2726ea7c1aeSJames Wright       }
273aa0b7f76SJames Wright     }
2744fa1625aSJames Wright   }
275aa0b7f76SJames Wright 
276aa0b7f76SJames Wright   // -- Compute and store anisotropy tensor
277e3663b90SJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_train_setup_data->elem_restr_grid_aniso,
278aa0b7f76SJames Wright                                                      &sgs_dd_train_setup_data->grid_aniso_ceed));
279aa0b7f76SJames Wright 
280aa0b7f76SJames Wright   // -- Create Nodal Evaluation Operator
281e3663b90SJames Wright   PetscCall(SetupTrainingDataCalculation(ceed, honee, problem, sgs_dd_train_setup_data));
282aa0b7f76SJames Wright 
283aa0b7f76SJames Wright   PetscCall(SGS_DD_TrainingSetupDataDestroy(sgs_dd_train_setup_data));
284aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
285aa0b7f76SJames Wright }
286aa0b7f76SJames Wright 
287aa0b7f76SJames Wright PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx) {
2880c373b74SJames Wright   Honee               honee = (Honee)ctx;
2890c373b74SJames Wright   Ceed                ceed  = honee->ceed;
29039169b57SJames Wright   SGS_DD_TrainingData sgs_dd_train;
2910c373b74SJames Wright   SmartSimData        smartsim = honee->smartsim;
292aa0b7f76SJames Wright   Vec                 TrainingData;
293ee81c423SRiccardo Balin   PetscMPIInt         rank;
294aa0b7f76SJames Wright 
295aa0b7f76SJames Wright   PetscFunctionBeginUser;
29639169b57SJames Wright   PetscCall(PetscObjectContainerQuery((PetscObject)honee, SGS_DD_TRAIN_KEY, &sgs_dd_train));
2970c373b74SJames Wright   PetscCallMPI(MPI_Comm_rank(honee->comm, &rank));
298ee81c423SRiccardo Balin 
299aa0b7f76SJames Wright   if (step_num % sgs_dd_train->write_data_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
300aa0b7f76SJames Wright   PetscCall(DMGetGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData));
301aa0b7f76SJames Wright 
3026ea7c1aeSJames Wright   for (PetscInt filter_index = 0; filter_index < sgs_dd_train->num_filter_widths; filter_index++) {
303ea615d4cSJames Wright     PetscCall(PetscLogEventBegin(HONEE_TrainDataCompute, 0, 0, 0, 0));
304aa0b7f76SJames Wright     {  // -- Compute and assemble training data
305aa0b7f76SJames Wright       Vec          FilteredVelocityGradient, FilteredFields, FilteredFields_loc;
306aa0b7f76SJames Wright       PetscMemType filtered_fields_mem_type;
307aa0b7f76SJames Wright       CeedVector   filtered_fields;
308aa0b7f76SJames Wright 
3096ea7c1aeSJames Wright       {  // Set filter width for the current solve
3106ea7c1aeSJames Wright         double       filter_width_scaling[3];
3116ea7c1aeSJames Wright         CeedOperator op_mat;
3126ea7c1aeSJames Wright         Mat          A_mat;
3136ea7c1aeSJames Wright 
3146ea7c1aeSJames Wright         for (int j = 0; j < 3; j++) filter_width_scaling[j] = sgs_dd_train->filter_widths[filter_index];
315*cb8a476cSJames Wright         PetscCall(KSPGetOperators(sgs_dd_train->diff_filter->ksp, &A_mat, NULL));
3166ea7c1aeSJames Wright         PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL));
317*cb8a476cSJames Wright         PetscCall(CeedOperatorSetContextDouble(op_mat, sgs_dd_train->diff_filter->filter_width_scaling_label, filter_width_scaling));
3186ea7c1aeSJames Wright       }
3196ea7c1aeSJames Wright 
320*cb8a476cSJames Wright       PetscCall(DMGetGlobalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields));
321*cb8a476cSJames Wright       PetscCall(DMGetLocalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields_loc));
322aa0b7f76SJames Wright 
323*cb8a476cSJames Wright       PetscCall(DifferentialFilterApply(honee, sgs_dd_train->diff_filter, solution_time, Q, FilteredFields));
324*cb8a476cSJames Wright       PetscCall(DMGlobalToLocal(sgs_dd_train->diff_filter->dm_filter, FilteredFields, INSERT_VALUES, FilteredFields_loc));
325aa0b7f76SJames Wright 
326aa0b7f76SJames Wright       PetscCall(DMGetGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient));
327aa0b7f76SJames Wright       PetscCall(VelocityGradientProjectionApply(sgs_dd_train->filtered_grad_velo_proj, FilteredFields_loc, FilteredVelocityGradient));
328aa0b7f76SJames Wright 
329aa0b7f76SJames Wright       {
330aa0b7f76SJames Wright         CeedOperatorField op_field;
3316ea7c1aeSJames Wright 
332aa0b7f76SJames Wright         PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->op_training_data_calc_ctx->op, "q", &op_field));
333aa0b7f76SJames Wright         PetscCallCeed(ceed, CeedOperatorFieldGetVector(op_field, &filtered_fields));
334aa0b7f76SJames Wright       }
3356ea7c1aeSJames Wright 
336a7dac1d5SJames Wright       PetscCall(VecPetscToCeed(FilteredFields_loc, &filtered_fields_mem_type, filtered_fields));  // filtered_fields is an implicit input
337aa0b7f76SJames Wright       PetscCall(ApplyCeedOperatorGlobalToGlobal(FilteredVelocityGradient, TrainingData, sgs_dd_train->op_training_data_calc_ctx));
338a7dac1d5SJames Wright       PetscCall(VecCeedToPetsc(filtered_fields, filtered_fields_mem_type, FilteredFields_loc));
339aa0b7f76SJames Wright 
340aa0b7f76SJames Wright       PetscCall(DMRestoreGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient));
341*cb8a476cSJames Wright       PetscCall(DMRestoreGlobalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields));
342*cb8a476cSJames Wright       PetscCall(DMRestoreLocalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields_loc));
343fff85bd3SJames Wright       PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
344aa0b7f76SJames Wright     }
345ea615d4cSJames Wright     PetscCall(PetscLogEventEnd(HONEE_TrainDataCompute, 0, 0, 0, 0));
346aa0b7f76SJames Wright 
347aa0b7f76SJames Wright     {  // -- Send training data to SmartSim
348aa0b7f76SJames Wright       char   array_key[PETSC_MAX_PATH_LEN];
349aa0b7f76SJames Wright       size_t array_key_len;
350aa0b7f76SJames Wright 
351aa0b7f76SJames Wright       if (sgs_dd_train->overwrite_training_data) {
3526ea7c1aeSJames Wright         PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index));
353aa0b7f76SJames Wright       } else {
3546ea7c1aeSJames Wright         PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index));
355aa0b7f76SJames Wright       }
356aa0b7f76SJames Wright       PetscCall(PetscStrlen(array_key, &array_key_len));
357aa0b7f76SJames Wright 
358aa0b7f76SJames Wright       {
359aa0b7f76SJames Wright         const PetscScalar *training_data;
360aa0b7f76SJames Wright         PetscCall(VecGetArrayRead(TrainingData, &training_data));
361ea615d4cSJames Wright         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Train, 0, 0, 0, 0));
36243e9749fSJames Wright         PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2,
363aa0b7f76SJames Wright                                        SRTensorTypeDouble, SRMemLayoutContiguous));
364ea615d4cSJames Wright         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Train, 0, 0, 0, 0));
365aa0b7f76SJames Wright         PetscCall(VecRestoreArrayRead(TrainingData, &training_data));
366aa0b7f76SJames Wright       }
367ee81c423SRiccardo Balin     }
368ee81c423SRiccardo Balin   }
369aa0b7f76SJames Wright 
370aa0b7f76SJames Wright   if (rank % smartsim->collocated_database_num_ranks == 0) {
3716ea7c1aeSJames Wright     const char tensor_name[] = "step";
372aa0b7f76SJames Wright     size_t     dim_2[1]      = {2};
373f4632befSRiccardo Balin     PetscInt64 step_array[2] = {step_num, step_num};
3746ea7c1aeSJames Wright 
375ea615d4cSJames Wright     PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
37643e9749fSJames Wright     PetscCallSmartRedis(
3776ea7c1aeSJames Wright         put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
378ea615d4cSJames Wright     PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
379aa0b7f76SJames Wright   }
380aa0b7f76SJames Wright 
38139169b57SJames Wright   PetscCall(DMRestoreGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData));
382aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
383aa0b7f76SJames Wright }
384aa0b7f76SJames Wright 
385632a41e1SJames Wright PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) {
3860c373b74SJames Wright   Honee        honee;
387632a41e1SJames Wright   const char   check_run_key[]   = "check-run";
388632a41e1SJames Wright   PetscReal    check_run[2]      = {1};
389632a41e1SJames Wright   const size_t check_run_dims[1] = {2};
390632a41e1SJames Wright   size_t       check_run_key_size;
391632a41e1SJames Wright 
392632a41e1SJames Wright   PetscFunctionBeginUser;
393632a41e1SJames Wright   PetscCall(PetscStrlen(check_run_key, &check_run_key_size));
3940c373b74SJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
3950c373b74SJames Wright   SmartSimData smartsim = honee->smartsim;
396632a41e1SJames Wright 
397ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
39843e9749fSJames Wright   PetscCallSmartRedis(
399632a41e1SJames Wright       unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous));
400ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
401632a41e1SJames Wright   if (check_run[0] == 0) {
4020c373b74SJames Wright     PetscCall(PetscPrintf(honee->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n"));
403632a41e1SJames Wright     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
404632a41e1SJames Wright   }
405632a41e1SJames Wright 
406632a41e1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
407632a41e1SJames Wright }
408