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
68d78d7c8SJames Wright #include <differential_filter.h>
7149fb536SJames Wright #include <navierstokes.h>
88d78d7c8SJames Wright #include <petscdmplex.h>
99ae013d6SJames Wright #include <smartsim-impl.h>
10aa0b7f76SJames Wright
11f5dc303cSJames Wright typedef struct SGS_DD_TrainingData_ *SGS_DD_TrainingData;
12f5dc303cSJames Wright struct SGS_DD_TrainingData_ {
1339169b57SJames Wright DM dm_dd_training;
1439169b57SJames Wright PetscInt num_comp_dd_inputs, write_data_interval, num_filter_widths;
1539169b57SJames Wright PetscScalar filter_widths[16];
1639169b57SJames Wright OperatorApplyContext op_training_data_calc_ctx;
17cb8a476cSJames Wright DiffFilterData diff_filter;
1839169b57SJames Wright NodalProjectionData filtered_grad_velo_proj;
1939169b57SJames Wright size_t training_data_array_dims[2];
2039169b57SJames Wright PetscBool overwrite_training_data;
21f5dc303cSJames Wright };
2239169b57SJames Wright
2339169b57SJames Wright #define SGS_DD_TRAIN_KEY "SGS Data Driven Training"
2439169b57SJames Wright
SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData * sgs_dd_train)259ae013d6SJames Wright static PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData *sgs_dd_train) {
2639169b57SJames Wright SGS_DD_TrainingData sgs_dd_train_ = *sgs_dd_train;
2739169b57SJames Wright
2839169b57SJames Wright PetscFunctionBeginUser;
2939169b57SJames Wright if (!sgs_dd_train_) PetscFunctionReturn(PETSC_SUCCESS);
3039169b57SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_train_->op_training_data_calc_ctx));
3139169b57SJames Wright PetscCall(NodalProjectionDataDestroy(&sgs_dd_train_->filtered_grad_velo_proj));
3239169b57SJames Wright PetscCall(DMDestroy(&sgs_dd_train_->dm_dd_training));
33cb8a476cSJames Wright PetscCall(DifferentialFilterDataDestroy(&sgs_dd_train_->diff_filter));
345206a5a0SJames Wright PetscCall(PetscFree(sgs_dd_train_));
355206a5a0SJames Wright *sgs_dd_train = NULL;
3639169b57SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
3739169b57SJames Wright }
3839169b57SJames Wright
3939169b57SJames Wright typedef struct {
40aa0b7f76SJames Wright CeedElemRestriction elem_restr_grid_aniso;
41aa0b7f76SJames Wright CeedVector grid_aniso_ceed;
42aa0b7f76SJames Wright CeedQFunctionContext sgs_dd_train_qfctx;
43aa0b7f76SJames Wright } *SGS_DD_TrainingSetupData;
44aa0b7f76SJames Wright
SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data)45aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
46aa0b7f76SJames Wright Ceed ceed;
47aa0b7f76SJames Wright
48aa0b7f76SJames Wright PetscFunctionBeginUser;
49aa0b7f76SJames Wright PetscCall(CeedElemRestrictionGetCeed(sgs_dd_train_setup_data->elem_restr_grid_aniso, &ceed));
50aa0b7f76SJames Wright
51aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_train_setup_data->elem_restr_grid_aniso));
52aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_train_setup_data->grid_aniso_ceed));
53aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_train_setup_data->sgs_dd_train_qfctx));
54aa0b7f76SJames Wright PetscCall(PetscFree(sgs_dd_train_setup_data));
55519781aeSJames Wright PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed");
56aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
57aa0b7f76SJames Wright }
58aa0b7f76SJames Wright
59aa0b7f76SJames Wright // @brief Create DM for storing data-drive SGS model inputs
SGS_DD_TrainingCreateDM(DM dm_source,DM * dm_dd_training,PetscInt degree,PetscInt q_extra,PetscInt * num_components)60aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
61aa0b7f76SJames Wright PetscSection section;
62aa0b7f76SJames Wright
63aa0b7f76SJames Wright PetscFunctionBeginUser;
64aa0b7f76SJames Wright *num_components = 12;
65aa0b7f76SJames Wright
66aa0b7f76SJames Wright PetscCall(DMClone(dm_source, dm_dd_training));
670dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(*dm_dd_training, PETSC_TRUE));
68aa0b7f76SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data"));
69aa0b7f76SJames Wright
70aa0b7f76SJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training));
71aa0b7f76SJames Wright
72aa0b7f76SJames Wright PetscCall(DMGetLocalSection(*dm_dd_training, §ion));
73aa0b7f76SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data"));
74aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1"));
75aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2"));
76aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3"));
77aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4"));
78aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5"));
79aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6"));
80aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX"));
81aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY"));
82aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ"));
83aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ"));
84aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ"));
85aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY"));
86aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
87aa0b7f76SJames Wright };
88aa0b7f76SJames Wright
89aa0b7f76SJames Wright // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes
SetupTrainingDataCalculation(Ceed ceed,Honee honee,ProblemData problem,SGS_DD_TrainingSetupData sgs_dd_train_setup_data)90e3663b90SJames Wright static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, Honee honee, ProblemData problem, SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
9139169b57SJames Wright SGS_DD_TrainingData sgs_dd_train;
92aa0b7f76SJames Wright CeedQFunction qf_sgs_dd_train;
93aa0b7f76SJames Wright CeedOperator op_sgs_dd_train;
94aa0b7f76SJames Wright CeedInt num_comp_grad_velo, num_comp_grid_aniso;
95aa0b7f76SJames Wright CeedVector inv_multiplicity, filtered_fields;
96aa0b7f76SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train;
97*e3db12f8SJames Wright PetscInt height = 0, dm_field = 0;
98aa0b7f76SJames Wright
99aa0b7f76SJames Wright PetscFunctionBeginUser;
1000c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_TRAIN_KEY, &sgs_dd_train));
101aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
102aa0b7f76SJames Wright
103*e3db12f8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field,
104*e3db12f8SJames Wright &elem_restr_sgs_train));
105*e3db12f8SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, PETSC_TRUE,
1065930f037SJames Wright &elem_restr_inv_multiplicity, &inv_multiplicity));
107aa0b7f76SJames Wright
108aa0b7f76SJames Wright CeedElemRestriction elem_restr_filtered_state;
109aa0b7f76SJames Wright CeedInt num_comp_filtered_state;
110aa0b7f76SJames Wright { // -- Setup filtered velocity gradient projection
111aa0b7f76SJames Wright CeedBasis basis_filtered_state;
112aa0b7f76SJames Wright CeedOperatorField op_field;
113cb8a476cSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->diff_filter->op_rhs_ctx->op, "v0", &op_field));
11401e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_filtered_state, &basis_filtered_state, NULL));
115aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state));
116e3663b90SJames Wright PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state,
117aa0b7f76SJames Wright &sgs_dd_train->filtered_grad_velo_proj));
118fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_filtered_state));
119aa0b7f76SJames Wright // Get velocity gradient information
120aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
121aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
122aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
123aa0b7f76SJames Wright }
124aa0b7f76SJames Wright
125aa0b7f76SJames Wright CeedElemRestriction elem_restr_filtered_velo_prod;
126aa0b7f76SJames Wright CeedInt num_comp_filtered_velo_prod;
127aa0b7f76SJames Wright { // Get filtered velocity product information
128aa0b7f76SJames Wright CeedOperatorField op_field;
129cb8a476cSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->diff_filter->op_rhs_ctx->op, "v1", &op_field));
130aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod));
131aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod));
132aa0b7f76SJames Wright }
133aa0b7f76SJames Wright
134aa0b7f76SJames Wright // -- Create operator for generating training data at nodes
135aa0b7f76SJames Wright // Differential Filter only provides filtered primitive variables
136aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim,
137aa0b7f76SJames Wright ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train));
138aa0b7f76SJames Wright
139aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx));
140aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE));
141aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE));
142aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
143aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
144aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE));
145aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE));
146aa0b7f76SJames Wright
147aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL));
148aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train));
149c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields));
150c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields));
151c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
152c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
153c6271fa9SJeremy L Thompson sgs_dd_train_setup_data->grid_aniso_ceed));
154c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
155c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
156aa0b7f76SJames Wright
157aa0b7f76SJames Wright PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL,
158aa0b7f76SJames Wright NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx));
159aa0b7f76SJames Wright
160aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
161aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
162aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
163fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_state));
164fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
165fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_velo_prod));
166aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_train));
167aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_train));
168aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
169aa0b7f76SJames Wright }
170aa0b7f76SJames Wright
SGS_DD_TrainingSetup(Ceed ceed,Honee honee)1719ae013d6SJames Wright PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee) {
1727f3a2123SJames Wright SGS_DDTrainingContext sgsdd_train_ctx;
173aa0b7f76SJames Wright SGS_DD_TrainingSetupData sgs_dd_train_setup_data;
17439169b57SJames Wright SGS_DD_TrainingData sgs_dd_train;
1759ae013d6SJames Wright ProblemData problem = honee->problem_data;
176aa0b7f76SJames Wright
177aa0b7f76SJames Wright PetscFunctionBeginUser;
1787f3a2123SJames Wright PetscCall(PetscNew(&sgsdd_train_ctx));
179aa0b7f76SJames Wright PetscCall(PetscNew(&sgs_dd_train_setup_data));
18039169b57SJames Wright PetscCall(PetscNew(&sgs_dd_train));
181f5dc303cSJames Wright *sgs_dd_train = (struct SGS_DD_TrainingData_){
182f5dc303cSJames Wright .overwrite_training_data = PETSC_TRUE,
183f5dc303cSJames Wright .write_data_interval = 1,
184f5dc303cSJames Wright .num_filter_widths = sizeof(sgs_dd_train->filter_widths) / sizeof(sgs_dd_train->filter_widths[0]),
185f5dc303cSJames Wright };
1860c70a8bcSJames Wright PetscCall(HoneeSetContainer(honee, SGS_DD_TRAIN_KEY, sgs_dd_train, (PetscCtxDestroyFn *)SGS_DD_TrainingDataDestroy));
187cb8a476cSJames Wright PetscCall(DifferentialFilterSetup(honee, &sgs_dd_train->diff_filter));
188aa0b7f76SJames Wright
1890c373b74SJames Wright PetscOptionsBegin(honee->comm, NULL, "SGS Data-Driven Training Options", NULL);
190aa0b7f76SJames Wright PetscCall(PetscOptionsInt("-sgs_train_write_data_interval", "Number of timesteps between writing data into database", NULL,
191aa0b7f76SJames Wright sgs_dd_train->write_data_interval, &sgs_dd_train->write_data_interval, NULL));
192aa0b7f76SJames Wright PetscCall(PetscOptionsBool("-sgs_train_overwrite_data", "Overwrite old training data in the database", NULL, sgs_dd_train->overwrite_training_data,
193aa0b7f76SJames Wright &sgs_dd_train->overwrite_training_data, NULL));
1946ea7c1aeSJames Wright PetscCall(PetscOptionsRealArray("-sgs_train_filter_width_scales", "Scales of each filter width put into training database", NULL,
1956ea7c1aeSJames Wright sgs_dd_train->filter_widths, &sgs_dd_train->num_filter_widths, NULL));
196aa0b7f76SJames Wright PetscOptionsEnd();
197aa0b7f76SJames Wright
198aa0b7f76SJames Wright // -- Create DM for storing training data
1990c373b74SJames Wright PetscCall(SGS_DD_TrainingCreateDM(honee->dm, &sgs_dd_train->dm_dd_training, honee->app_ctx->degree, honee->app_ctx->q_extra,
200aa0b7f76SJames Wright &sgs_dd_train->num_comp_dd_inputs));
201aa0b7f76SJames Wright
202aa0b7f76SJames Wright { // -- Create QFunction Context
203cde3d787SJames Wright NewtonianIdealGasContext newt_ctx;
204cde3d787SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &newt_ctx));
205cde3d787SJames Wright sgsdd_train_ctx->newt_ctx = *newt_ctx;
206cde3d787SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &newt_ctx));
2070c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_train_setup_data->sgs_dd_train_qfctx));
208aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, CEED_USE_POINTER,
2097f3a2123SJames Wright sizeof(*sgsdd_train_ctx), sgsdd_train_ctx));
210aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, FreeContextPetsc));
211aa0b7f76SJames Wright }
212aa0b7f76SJames Wright
213aa0b7f76SJames Wright { // -- Send training data array info to SmartRedis database
214aa0b7f76SJames Wright PetscMPIInt rank, num_ranks;
2157ebeccb9SJames Wright SmartSimData smartsim;
2167ebeccb9SJames Wright PetscCall(HoneeGetSmartSimData(honee, &smartsim));
2170c373b74SJames Wright PetscCallMPI(MPI_Comm_rank(honee->comm, &rank));
2180c373b74SJames Wright PetscCallMPI(MPI_Comm_size(honee->comm, &num_ranks));
219aa0b7f76SJames Wright
220aa0b7f76SJames Wright {
221aa0b7f76SJames Wright PetscSection global_section;
2227ff16c02SJames Wright PetscInt num_dofs, num_comps, local_min_max[2] = {0.}, global_min_max[2] = {0.};
2237ff16c02SJames Wright
224aa0b7f76SJames Wright PetscCall(DMGetGlobalSection(sgs_dd_train->dm_dd_training, &global_section));
225aa0b7f76SJames Wright PetscCall(DMGetGlobalVectorInfo(sgs_dd_train->dm_dd_training, &num_dofs, NULL, NULL));
226aa0b7f76SJames Wright PetscCall(PetscSectionGetFieldComponents(global_section, 0, &num_comps));
2277ff16c02SJames Wright local_min_max[0] = num_dofs;
2280c373b74SJames Wright PetscCall(PetscGlobalMinMaxInt(honee->comm, local_min_max, global_min_max));
2297ff16c02SJames Wright
2307ff16c02SJames Wright sgs_dd_train->training_data_array_dims[0] = global_min_max[0] / num_comps;
231aa0b7f76SJames Wright sgs_dd_train->training_data_array_dims[1] = num_comps;
232aa0b7f76SJames Wright }
233aa0b7f76SJames Wright
234aa0b7f76SJames Wright if (rank % smartsim->collocated_database_num_ranks == 0) {
2354fa1625aSJames Wright { // Communicate info on simulation size
2364fa1625aSJames Wright const char tensor_name[] = "sizeInfo";
237aa0b7f76SJames Wright size_t array_info_dim = 6;
238f4632befSRiccardo Balin PetscInt64 array_info[6] = {0}, num_features = 6;
239aa0b7f76SJames Wright
240aa0b7f76SJames Wright array_info[0] = sgs_dd_train->training_data_array_dims[0];
241aa0b7f76SJames Wright array_info[1] = sgs_dd_train->training_data_array_dims[1];
242aa0b7f76SJames Wright array_info[2] = num_features;
243aa0b7f76SJames Wright array_info[3] = num_ranks;
244aa0b7f76SJames Wright array_info[4] = smartsim->collocated_database_num_ranks;
245aa0b7f76SJames Wright array_info[5] = rank;
246aa0b7f76SJames Wright
247ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
2489eadbee4SJames Wright PetscCallSmartRedis(put_tensor(smartsim->client, tensor_name, strlen(tensor_name), array_info, &array_info_dim, 1, SRTensorTypeInt64,
2499eadbee4SJames Wright SRMemLayoutContiguous));
2504fa1625aSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
251ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
2524fa1625aSJames Wright }
253aa0b7f76SJames Wright
2544fa1625aSJames Wright { // Send array that communicates if tensors are overwritten in database
2554fa1625aSJames Wright const char tensor_name[] = "tensor-ow";
256f4632befSRiccardo Balin PetscInt64 tensor_overwrite[2] = {sgs_dd_train->overwrite_training_data};
257aa0b7f76SJames Wright size_t dim_2[1] = {2};
2584fa1625aSJames Wright
259ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
2609eadbee4SJames Wright PetscCallSmartRedis(put_tensor(smartsim->client, tensor_name, strlen(tensor_name), tensor_overwrite, dim_2, 1, SRTensorTypeInt64,
2619eadbee4SJames Wright SRMemLayoutContiguous));
2624fa1625aSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
263ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
264aa0b7f76SJames Wright }
2656ea7c1aeSJames Wright
2666ea7c1aeSJames Wright { // Communicate number of filter widths used
2676ea7c1aeSJames Wright const char tensor_name[] = "num_filter_widths";
2686ea7c1aeSJames Wright PetscInt64 num_filter_widths = sgs_dd_train->num_filter_widths;
2696ea7c1aeSJames Wright size_t dim_2 = 1;
2706ea7c1aeSJames Wright
271ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
2729eadbee4SJames Wright PetscCallSmartRedis(put_tensor(smartsim->client, tensor_name, strlen(tensor_name), &num_filter_widths, &dim_2, 1, SRTensorTypeInt64,
2739eadbee4SJames Wright SRMemLayoutContiguous));
2746ea7c1aeSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
275ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
2766ea7c1aeSJames Wright }
277aa0b7f76SJames Wright }
2784fa1625aSJames Wright }
279aa0b7f76SJames Wright
280aa0b7f76SJames Wright // -- Compute and store anisotropy tensor
281e3663b90SJames Wright PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_train_setup_data->elem_restr_grid_aniso,
282aa0b7f76SJames Wright &sgs_dd_train_setup_data->grid_aniso_ceed));
283aa0b7f76SJames Wright
284aa0b7f76SJames Wright // -- Create Nodal Evaluation Operator
285e3663b90SJames Wright PetscCall(SetupTrainingDataCalculation(ceed, honee, problem, sgs_dd_train_setup_data));
286aa0b7f76SJames Wright
287aa0b7f76SJames Wright PetscCall(SGS_DD_TrainingSetupDataDestroy(sgs_dd_train_setup_data));
288aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
289aa0b7f76SJames Wright }
290aa0b7f76SJames Wright
TSMonitor_SGS_DD_Training(TS ts,PetscInt step_num,PetscReal solution_time,Vec Q,void * ctx)291aa0b7f76SJames Wright PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx) {
2920c373b74SJames Wright Honee honee = (Honee)ctx;
2930c373b74SJames Wright Ceed ceed = honee->ceed;
29439169b57SJames Wright SGS_DD_TrainingData sgs_dd_train;
2957ebeccb9SJames Wright SmartSimData smartsim;
296aa0b7f76SJames Wright Vec TrainingData;
297ee81c423SRiccardo Balin PetscMPIInt rank;
298aa0b7f76SJames Wright
299aa0b7f76SJames Wright PetscFunctionBeginUser;
3007ebeccb9SJames Wright PetscCall(HoneeGetSmartSimData(honee, &smartsim));
3010c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_TRAIN_KEY, &sgs_dd_train));
3020c373b74SJames Wright PetscCallMPI(MPI_Comm_rank(honee->comm, &rank));
303ee81c423SRiccardo Balin
304aa0b7f76SJames Wright if (step_num % sgs_dd_train->write_data_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
305aa0b7f76SJames Wright PetscCall(DMGetGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData));
306aa0b7f76SJames Wright
3076ea7c1aeSJames Wright for (PetscInt filter_index = 0; filter_index < sgs_dd_train->num_filter_widths; filter_index++) {
308ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_TrainDataCompute, 0, 0, 0, 0));
309aa0b7f76SJames Wright { // -- Compute and assemble training data
310aa0b7f76SJames Wright Vec FilteredVelocityGradient, FilteredFields, FilteredFields_loc;
311aa0b7f76SJames Wright PetscMemType filtered_fields_mem_type;
312aa0b7f76SJames Wright CeedVector filtered_fields;
313aa0b7f76SJames Wright
3146ea7c1aeSJames Wright { // Set filter width for the current solve
3156ea7c1aeSJames Wright double filter_width_scaling[3];
3166ea7c1aeSJames Wright CeedOperator op_mat;
3176ea7c1aeSJames Wright Mat A_mat;
3186ea7c1aeSJames Wright
3196ea7c1aeSJames Wright for (int j = 0; j < 3; j++) filter_width_scaling[j] = sgs_dd_train->filter_widths[filter_index];
320cb8a476cSJames Wright PetscCall(KSPGetOperators(sgs_dd_train->diff_filter->ksp, &A_mat, NULL));
3216ea7c1aeSJames Wright PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL));
322cb8a476cSJames Wright PetscCall(CeedOperatorSetContextDouble(op_mat, sgs_dd_train->diff_filter->filter_width_scaling_label, filter_width_scaling));
3236ea7c1aeSJames Wright }
3246ea7c1aeSJames Wright
325cb8a476cSJames Wright PetscCall(DMGetGlobalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields));
326cb8a476cSJames Wright PetscCall(DMGetLocalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields_loc));
327aa0b7f76SJames Wright
328cb8a476cSJames Wright PetscCall(DifferentialFilterApply(honee, sgs_dd_train->diff_filter, solution_time, Q, FilteredFields));
329cb8a476cSJames Wright PetscCall(DMGlobalToLocal(sgs_dd_train->diff_filter->dm_filter, FilteredFields, INSERT_VALUES, FilteredFields_loc));
330aa0b7f76SJames Wright
331aa0b7f76SJames Wright PetscCall(DMGetGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient));
332aa0b7f76SJames Wright PetscCall(VelocityGradientProjectionApply(sgs_dd_train->filtered_grad_velo_proj, FilteredFields_loc, FilteredVelocityGradient));
333aa0b7f76SJames Wright
334aa0b7f76SJames Wright {
335aa0b7f76SJames Wright CeedOperatorField op_field;
3366ea7c1aeSJames Wright
337aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->op_training_data_calc_ctx->op, "q", &op_field));
338aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(op_field, &filtered_fields));
339aa0b7f76SJames Wright }
3406ea7c1aeSJames Wright
341a7dac1d5SJames Wright PetscCall(VecPetscToCeed(FilteredFields_loc, &filtered_fields_mem_type, filtered_fields)); // filtered_fields is an implicit input
342aa0b7f76SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(FilteredVelocityGradient, TrainingData, sgs_dd_train->op_training_data_calc_ctx));
343a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(filtered_fields, filtered_fields_mem_type, FilteredFields_loc));
344aa0b7f76SJames Wright
345aa0b7f76SJames Wright PetscCall(DMRestoreGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient));
346cb8a476cSJames Wright PetscCall(DMRestoreGlobalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields));
347cb8a476cSJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields_loc));
348fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
349aa0b7f76SJames Wright }
350ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_TrainDataCompute, 0, 0, 0, 0));
351aa0b7f76SJames Wright
352aa0b7f76SJames Wright { // -- Send training data to SmartSim
353aa0b7f76SJames Wright char array_key[PETSC_MAX_PATH_LEN];
354aa0b7f76SJames Wright size_t array_key_len;
355aa0b7f76SJames Wright
356aa0b7f76SJames Wright if (sgs_dd_train->overwrite_training_data) {
3576ea7c1aeSJames Wright PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index));
358aa0b7f76SJames Wright } else {
3596ea7c1aeSJames Wright PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index));
360aa0b7f76SJames Wright }
361aa0b7f76SJames Wright PetscCall(PetscStrlen(array_key, &array_key_len));
362aa0b7f76SJames Wright
363aa0b7f76SJames Wright {
364aa0b7f76SJames Wright const PetscScalar *training_data;
365aa0b7f76SJames Wright PetscCall(VecGetArrayRead(TrainingData, &training_data));
366e171caa6SJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Write, 0, 0, 0, 0));
36743e9749fSJames Wright PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2,
368aa0b7f76SJames Wright SRTensorTypeDouble, SRMemLayoutContiguous));
369e171caa6SJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Write, 0, 0, 0, 0));
370aa0b7f76SJames Wright PetscCall(VecRestoreArrayRead(TrainingData, &training_data));
371aa0b7f76SJames Wright }
372ee81c423SRiccardo Balin }
373ee81c423SRiccardo Balin }
374aa0b7f76SJames Wright
375aa0b7f76SJames Wright if (rank % smartsim->collocated_database_num_ranks == 0) {
3766ea7c1aeSJames Wright const char tensor_name[] = "step";
377aa0b7f76SJames Wright size_t dim_2[1] = {2};
378f4632befSRiccardo Balin PetscInt64 step_array[2] = {step_num, step_num};
3796ea7c1aeSJames Wright
380ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
3819eadbee4SJames Wright PetscCallSmartRedis(put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64,
3829eadbee4SJames Wright SRMemLayoutContiguous));
383ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
384aa0b7f76SJames Wright }
385aa0b7f76SJames Wright
38639169b57SJames Wright PetscCall(DMRestoreGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData));
387aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
388aa0b7f76SJames Wright }
389aa0b7f76SJames Wright
TSPostStep_SGS_DD_Training(TS ts)390632a41e1SJames Wright PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) {
3910c373b74SJames Wright Honee honee;
392632a41e1SJames Wright const char check_run_key[] = "check-run";
393632a41e1SJames Wright PetscReal check_run[2] = {1};
394632a41e1SJames Wright const size_t check_run_dims[1] = {2};
395632a41e1SJames Wright size_t check_run_key_size;
3967ebeccb9SJames Wright SmartSimData smartsim;
397632a41e1SJames Wright
398632a41e1SJames Wright PetscFunctionBeginUser;
399632a41e1SJames Wright PetscCall(PetscStrlen(check_run_key, &check_run_key_size));
4000c373b74SJames Wright PetscCall(TSGetApplicationContext(ts, &honee));
4017ebeccb9SJames Wright PetscCall(HoneeGetSmartSimData(honee, &smartsim));
402632a41e1SJames Wright
403ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
4049eadbee4SJames Wright PetscCallSmartRedis(unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble,
4059eadbee4SJames Wright SRMemLayoutContiguous));
406ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
407632a41e1SJames Wright if (check_run[0] == 0) {
4080c373b74SJames Wright PetscCall(PetscPrintf(honee->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n"));
409632a41e1SJames Wright PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
410632a41e1SJames Wright }
411632a41e1SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
412632a41e1SJames Wright }
413