xref: /honee/src/smartsim/sgs_dd_training.c (revision ea615d4cc464aa6ad650c06fae6d120cc2465bc4)
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 {
12aa0b7f76SJames Wright   CeedElemRestriction  elem_restr_grid_aniso;
13aa0b7f76SJames Wright   CeedVector           grid_aniso_ceed;
14aa0b7f76SJames Wright   CeedQFunctionContext sgs_dd_train_qfctx;
15aa0b7f76SJames Wright } *SGS_DD_TrainingSetupData;
16aa0b7f76SJames Wright 
17aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
18aa0b7f76SJames Wright   Ceed ceed;
19aa0b7f76SJames Wright 
20aa0b7f76SJames Wright   PetscFunctionBeginUser;
21aa0b7f76SJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_train_setup_data->elem_restr_grid_aniso, &ceed));
22aa0b7f76SJames Wright 
23aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_train_setup_data->elem_restr_grid_aniso));
24aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_train_setup_data->grid_aniso_ceed));
25aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_train_setup_data->sgs_dd_train_qfctx));
26aa0b7f76SJames Wright   PetscCall(PetscFree(sgs_dd_train_setup_data));
27519781aeSJames Wright   PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed");
28aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
29aa0b7f76SJames Wright }
30aa0b7f76SJames Wright 
31aa0b7f76SJames Wright // @brief Create DM for storing data-drive SGS model inputs
32aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
33aa0b7f76SJames Wright   PetscSection section;
34aa0b7f76SJames Wright 
35aa0b7f76SJames Wright   PetscFunctionBeginUser;
36aa0b7f76SJames Wright   *num_components = 12;
37aa0b7f76SJames Wright 
38aa0b7f76SJames Wright   PetscCall(DMClone(dm_source, dm_dd_training));
390dee9b8eSJames Wright   PetscCall(DMSetMatrixPreallocateSkip(*dm_dd_training, PETSC_TRUE));
40aa0b7f76SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data"));
41aa0b7f76SJames Wright 
42aa0b7f76SJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training));
43aa0b7f76SJames Wright 
44aa0b7f76SJames Wright   PetscCall(DMGetLocalSection(*dm_dd_training, &section));
45aa0b7f76SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data"));
46aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1"));
47aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2"));
48aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3"));
49aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4"));
50aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5"));
51aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6"));
52aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX"));
53aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY"));
54aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ"));
55aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ"));
56aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ"));
57aa0b7f76SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY"));
58aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
59aa0b7f76SJames Wright };
60aa0b7f76SJames Wright 
61aa0b7f76SJames Wright // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes
62e3663b90SJames Wright static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, Honee honee, ProblemData problem, SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
630c373b74SJames Wright   SGS_DD_TrainingData sgs_dd_train = honee->sgs_dd_train;
64aa0b7f76SJames Wright   CeedQFunction       qf_sgs_dd_train;
65aa0b7f76SJames Wright   CeedOperator        op_sgs_dd_train;
66aa0b7f76SJames Wright   CeedInt             num_comp_grad_velo, num_comp_grid_aniso;
67aa0b7f76SJames Wright   CeedVector          inv_multiplicity, filtered_fields;
68aa0b7f76SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train;
69aa0b7f76SJames Wright   DMLabel             domain_label = NULL;
70aa0b7f76SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
71aa0b7f76SJames Wright 
72aa0b7f76SJames Wright   PetscFunctionBeginUser;
73aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
74aa0b7f76SJames Wright 
75aa0b7f76SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, &elem_restr_sgs_train));
765930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, PETSC_TRUE,
775930f037SJames Wright                                    &elem_restr_inv_multiplicity, &inv_multiplicity));
78aa0b7f76SJames Wright 
79aa0b7f76SJames Wright   CeedElemRestriction elem_restr_filtered_state;
80aa0b7f76SJames Wright   CeedInt             num_comp_filtered_state;
81aa0b7f76SJames Wright   {  // -- Setup filtered velocity gradient projection
82aa0b7f76SJames Wright     CeedBasis         basis_filtered_state;
83aa0b7f76SJames Wright     CeedOperatorField op_field;
840c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->diff_filter->op_rhs_ctx->op, "v0", &op_field));
8501e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_filtered_state, &basis_filtered_state, NULL));
86aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state));
87e3663b90SJames Wright     PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state,
88aa0b7f76SJames Wright                                               &sgs_dd_train->filtered_grad_velo_proj));
89fff85bd3SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_filtered_state));
90aa0b7f76SJames Wright     // Get velocity gradient information
91aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
92aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
93aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
94aa0b7f76SJames Wright   }
95aa0b7f76SJames Wright 
96aa0b7f76SJames Wright   CeedElemRestriction elem_restr_filtered_velo_prod;
97aa0b7f76SJames Wright   CeedInt             num_comp_filtered_velo_prod;
98aa0b7f76SJames Wright   {  // Get filtered velocity product information
99aa0b7f76SJames Wright     CeedOperatorField op_field;
1000c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->diff_filter->op_rhs_ctx->op, "v1", &op_field));
101aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod));
102aa0b7f76SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod));
103aa0b7f76SJames Wright   }
104aa0b7f76SJames Wright 
105aa0b7f76SJames Wright   // -- Create operator for generating training data at nodes
106aa0b7f76SJames Wright   // Differential Filter only provides filtered primitive variables
107aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim,
108aa0b7f76SJames Wright                                                   ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train));
109aa0b7f76SJames Wright 
110aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx));
111aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE));
112aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE));
113aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
114aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
115aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE));
116aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE));
117aa0b7f76SJames Wright 
118aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL));
119aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train));
120c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields));
121c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields));
122c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
123c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
124c6271fa9SJeremy L Thompson                                            sgs_dd_train_setup_data->grid_aniso_ceed));
125c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
126c6271fa9SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
127aa0b7f76SJames Wright 
128aa0b7f76SJames Wright   PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL,
129aa0b7f76SJames Wright                                        NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx));
130aa0b7f76SJames Wright 
131aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
132aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
133aa0b7f76SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
134fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_state));
135fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
136fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_velo_prod));
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 
142e3663b90SJames Wright PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem) {
143aa0b7f76SJames Wright   SGS_DDTrainingContext    sgsdd_train_qfctx;
144aa0b7f76SJames Wright   SGS_DD_TrainingSetupData sgs_dd_train_setup_data;
145aa0b7f76SJames Wright 
146aa0b7f76SJames Wright   PetscFunctionBeginUser;
147e3663b90SJames Wright   if (!honee->diff_filter) PetscCall(DifferentialFilterSetup(ceed, honee, problem));
1480c373b74SJames Wright   if (!honee->smartsim) PetscCall(SmartSimSetup(honee));
149aa0b7f76SJames Wright 
150aa0b7f76SJames Wright   PetscCall(PetscNew(&sgsdd_train_qfctx));
151aa0b7f76SJames Wright   PetscCall(PetscNew(&sgs_dd_train_setup_data));
1520c373b74SJames Wright   PetscCall(PetscNew(&honee->sgs_dd_train));
1530c373b74SJames Wright   SGS_DD_TrainingData sgs_dd_train = honee->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]);
1580c373b74SJames Wright   PetscOptionsBegin(honee->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
1680c373b74SJames Wright   PetscCall(SGS_DD_TrainingCreateDM(honee->dm, &sgs_dd_train->dm_dd_training, honee->app_ctx->degree, honee->app_ctx->q_extra,
169aa0b7f76SJames Wright                                     &sgs_dd_train->num_comp_dd_inputs));
170aa0b7f76SJames Wright 
171aa0b7f76SJames Wright   {  // -- Create QFunction Context
172aa0b7f76SJames Wright     NewtonianIdealGasContext gas;
173e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas));
174aa0b7f76SJames Wright     sgsdd_train_qfctx->gas = *gas;
175e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas));
1760c373b74SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->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;
1840c373b74SJames Wright     SmartSimData smartsim = honee->smartsim;
1850c373b74SJames Wright     PetscCallMPI(MPI_Comm_rank(honee->comm, &rank));
1860c373b74SJames Wright     PetscCallMPI(MPI_Comm_size(honee->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;
1960c373b74SJames Wright       PetscCall(PetscGlobalMinMaxInt(honee->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 
215*ea615d4cSJames Wright         PetscCall(PetscLogEventBegin(HONEE_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)));
219*ea615d4cSJames Wright         PetscCall(PetscLogEventEnd(HONEE_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 
227*ea615d4cSJames Wright         PetscCall(PetscLogEventBegin(HONEE_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)));
231*ea615d4cSJames Wright         PetscCall(PetscLogEventEnd(HONEE_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 
239*ea615d4cSJames Wright         PetscCall(PetscLogEventBegin(HONEE_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)));
243*ea615d4cSJames Wright         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
2446ea7c1aeSJames Wright       }
245aa0b7f76SJames Wright     }
2464fa1625aSJames Wright   }
247aa0b7f76SJames Wright 
248aa0b7f76SJames Wright   // -- Compute and store anisotropy tensor
249e3663b90SJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &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
253e3663b90SJames Wright   PetscCall(SetupTrainingDataCalculation(ceed, honee, 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) {
2600c373b74SJames Wright   Honee               honee        = (Honee)ctx;
2610c373b74SJames Wright   Ceed                ceed         = honee->ceed;
2620c373b74SJames Wright   SGS_DD_TrainingData sgs_dd_train = honee->sgs_dd_train;
2630c373b74SJames Wright   SmartSimData        smartsim     = honee->smartsim;
264aa0b7f76SJames Wright   Vec                 TrainingData;
265ee81c423SRiccardo Balin   PetscMPIInt         rank;
266aa0b7f76SJames Wright 
267aa0b7f76SJames Wright   PetscFunctionBeginUser;
268ee81c423SRiccardo Balin 
2690c373b74SJames Wright   PetscCallMPI(MPI_Comm_rank(honee->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++) {
275*ea615d4cSJames Wright     PetscCall(PetscLogEventBegin(HONEE_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];
2870c373b74SJames Wright         PetscCall(KSPGetOperators(honee->diff_filter->ksp, &A_mat, NULL));
2886ea7c1aeSJames Wright         PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL));
2890c373b74SJames Wright         PetscCall(CeedOperatorSetContextDouble(op_mat, honee->diff_filter->filter_width_scaling_label, filter_width_scaling));
2906ea7c1aeSJames Wright       }
2916ea7c1aeSJames Wright 
2920c373b74SJames Wright       PetscCall(DMGetGlobalVector(honee->diff_filter->dm_filter, &FilteredFields));
2930c373b74SJames Wright       PetscCall(DMGetLocalVector(honee->diff_filter->dm_filter, &FilteredFields_loc));
294aa0b7f76SJames Wright 
2950c373b74SJames Wright       PetscCall(DifferentialFilterApply(honee, solution_time, Q, FilteredFields));
2960c373b74SJames Wright       PetscCall(DMGlobalToLocal(honee->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));
3130c373b74SJames Wright       PetscCall(DMRestoreGlobalVector(honee->diff_filter->dm_filter, &FilteredFields));
3140c373b74SJames Wright       PetscCall(DMRestoreLocalVector(honee->diff_filter->dm_filter, &FilteredFields_loc));
315fff85bd3SJames Wright       PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
316aa0b7f76SJames Wright     }
317*ea615d4cSJames Wright     PetscCall(PetscLogEventEnd(HONEE_TrainDataCompute, 0, 0, 0, 0));
318aa0b7f76SJames Wright 
319aa0b7f76SJames Wright     {  // -- Send training data to SmartSim
320aa0b7f76SJames Wright       char   array_key[PETSC_MAX_PATH_LEN];
321aa0b7f76SJames Wright       size_t array_key_len;
322aa0b7f76SJames Wright 
323aa0b7f76SJames Wright       if (sgs_dd_train->overwrite_training_data) {
3246ea7c1aeSJames Wright         PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index));
325aa0b7f76SJames Wright       } else {
3266ea7c1aeSJames Wright         PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index));
327aa0b7f76SJames Wright       }
328aa0b7f76SJames Wright       PetscCall(PetscStrlen(array_key, &array_key_len));
329aa0b7f76SJames Wright 
330aa0b7f76SJames Wright       {
331aa0b7f76SJames Wright         const PetscScalar *training_data;
332aa0b7f76SJames Wright         PetscCall(VecGetArrayRead(TrainingData, &training_data));
333*ea615d4cSJames Wright         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Train, 0, 0, 0, 0));
33443e9749fSJames Wright         PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2,
335aa0b7f76SJames Wright                                        SRTensorTypeDouble, SRMemLayoutContiguous));
336*ea615d4cSJames Wright         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Train, 0, 0, 0, 0));
337aa0b7f76SJames Wright         PetscCall(VecRestoreArrayRead(TrainingData, &training_data));
338aa0b7f76SJames Wright       }
339ee81c423SRiccardo Balin     }
340ee81c423SRiccardo Balin   }
341aa0b7f76SJames Wright 
342aa0b7f76SJames Wright   if (rank % smartsim->collocated_database_num_ranks == 0) {
3436ea7c1aeSJames Wright     const char tensor_name[] = "step";
344aa0b7f76SJames Wright     size_t     dim_2[1]      = {2};
345f4632befSRiccardo Balin     PetscInt64 step_array[2] = {step_num, step_num};
3466ea7c1aeSJames Wright 
347*ea615d4cSJames Wright     PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
34843e9749fSJames Wright     PetscCallSmartRedis(
3496ea7c1aeSJames Wright         put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
350*ea615d4cSJames Wright     PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
351aa0b7f76SJames Wright   }
352aa0b7f76SJames Wright 
3530c373b74SJames Wright   PetscCall(DMRestoreGlobalVector(honee->sgs_dd_train->dm_dd_training, &TrainingData));
354aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
355aa0b7f76SJames Wright }
356aa0b7f76SJames Wright 
357632a41e1SJames Wright PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) {
3580c373b74SJames Wright   Honee        honee;
359632a41e1SJames Wright   const char   check_run_key[]   = "check-run";
360632a41e1SJames Wright   PetscReal    check_run[2]      = {1};
361632a41e1SJames Wright   const size_t check_run_dims[1] = {2};
362632a41e1SJames Wright   size_t       check_run_key_size;
363632a41e1SJames Wright 
364632a41e1SJames Wright   PetscFunctionBeginUser;
365632a41e1SJames Wright   PetscCall(PetscStrlen(check_run_key, &check_run_key_size));
3660c373b74SJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
3670c373b74SJames Wright   SmartSimData smartsim = honee->smartsim;
368632a41e1SJames Wright 
369*ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
37043e9749fSJames Wright   PetscCallSmartRedis(
371632a41e1SJames Wright       unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous));
372*ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
373632a41e1SJames Wright   if (check_run[0] == 0) {
3740c373b74SJames Wright     PetscCall(PetscPrintf(honee->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n"));
375632a41e1SJames Wright     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
376632a41e1SJames Wright   }
377632a41e1SJames Wright 
378632a41e1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
379632a41e1SJames Wright }
380632a41e1SJames Wright 
381aa0b7f76SJames Wright PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train) {
382aa0b7f76SJames Wright   PetscFunctionBeginUser;
383aa0b7f76SJames Wright   if (!sgs_dd_train) PetscFunctionReturn(PETSC_SUCCESS);
384aa0b7f76SJames Wright 
385aa0b7f76SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_train->op_training_data_calc_ctx));
386aa0b7f76SJames Wright   PetscCall(NodalProjectionDataDestroy(sgs_dd_train->filtered_grad_velo_proj));
387aa0b7f76SJames Wright   PetscCall(DMDestroy(&sgs_dd_train->dm_dd_training));
388aa0b7f76SJames Wright   PetscCall(PetscFree(sgs_dd_train));
389aa0b7f76SJames Wright 
390aa0b7f76SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
391aa0b7f76SJames Wright }
392