xref: /honee/src/smartsim/solution.c (revision 9ae013d64647ae6ed3b5a636462cb64250294c78)
1*9ae013d6SJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2025, HONEE contributors.
2*9ae013d6SJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3*9ae013d6SJames Wright 
4*9ae013d6SJames Wright #include <navierstokes.h>
5*9ae013d6SJames Wright #include <smartsim-impl.h>
6*9ae013d6SJames Wright 
7*9ae013d6SJames Wright typedef struct {
8*9ae013d6SJames Wright   PetscInt  write_data_interval;
9*9ae013d6SJames Wright   size_t    local_array_dims[2];
10*9ae013d6SJames Wright   PetscBool overwrite_training_data;
11*9ae013d6SJames Wright } *SmartSimSolutionData;
12*9ae013d6SJames Wright 
13*9ae013d6SJames Wright PetscErrorCode TSMonitor_SmartSimSolutionSetup(TS ts, PetscViewerAndFormat *ctx) {
14*9ae013d6SJames Wright   SmartSimSolutionData smartsimsol;
15*9ae013d6SJames Wright   Honee                honee;
16*9ae013d6SJames Wright   MPI_Comm             comm = PetscObjectComm((PetscObject)ts);
17*9ae013d6SJames Wright 
18*9ae013d6SJames Wright   PetscFunctionBeginUser;
19*9ae013d6SJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
20*9ae013d6SJames Wright   PetscCall(PetscNew(&smartsimsol));
21*9ae013d6SJames Wright   smartsimsol->overwrite_training_data = PETSC_TRUE;
22*9ae013d6SJames Wright 
23*9ae013d6SJames Wright   PetscOptionsBegin(comm, NULL, "SmartSim Solution Writing", NULL);
24*9ae013d6SJames Wright   PetscCall(PetscOptionsBool("-ts_monitor_smartsim_solution_overwrite_data", "Overwrite old solution data in the database", NULL,
25*9ae013d6SJames Wright                              smartsimsol->overwrite_training_data, &smartsimsol->overwrite_training_data, NULL));
26*9ae013d6SJames Wright   PetscOptionsEnd();
27*9ae013d6SJames Wright 
28*9ae013d6SJames Wright   {  // Get solution vector size
29*9ae013d6SJames Wright     PetscSection output_section;
30*9ae013d6SJames Wright     PetscInt     num_dofs, num_comps;
31*9ae013d6SJames Wright     DM           dm, output_dm;
32*9ae013d6SJames Wright 
33*9ae013d6SJames Wright     PetscCall(TSGetDM(ts, &dm));
34*9ae013d6SJames Wright     PetscCall(DMGetOutputDM(dm, &output_dm));
35*9ae013d6SJames Wright 
36*9ae013d6SJames Wright     PetscCall(DMGetGlobalVectorInfo(output_dm, &num_dofs, NULL, NULL));
37*9ae013d6SJames Wright     PetscCall(DMGetGlobalSection(output_dm, &output_section));
38*9ae013d6SJames Wright     PetscCall(PetscSectionGetFieldComponents(output_section, 0, &num_comps));
39*9ae013d6SJames Wright     smartsimsol->local_array_dims[0] = num_dofs / num_comps;
40*9ae013d6SJames Wright     smartsimsol->local_array_dims[1] = num_comps;
41*9ae013d6SJames Wright   }
42*9ae013d6SJames Wright   ctx->data = smartsimsol;
43*9ae013d6SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
44*9ae013d6SJames Wright }
45*9ae013d6SJames Wright 
46*9ae013d6SJames Wright PetscErrorCode TSMonitor_SmartSimSolution(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) {
47*9ae013d6SJames Wright   Honee        honee;
48*9ae013d6SJames Wright   SmartSimData smartsim;
49*9ae013d6SJames Wright   Vec          Q_output;
50*9ae013d6SJames Wright   PetscMPIInt  rank;
51*9ae013d6SJames Wright   DM           dm, output_dm;
52*9ae013d6SJames Wright 
53*9ae013d6SJames Wright   PetscFunctionBeginUser;
54*9ae013d6SJames Wright   SmartSimSolutionData smartsimsol = ctx->data;
55*9ae013d6SJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
56*9ae013d6SJames Wright   PetscCall(HoneeGetSmartSimData(honee, &smartsim));
57*9ae013d6SJames Wright   PetscCallMPI(MPI_Comm_rank(honee->comm, &rank));
58*9ae013d6SJames Wright 
59*9ae013d6SJames Wright   if (step_num % ctx->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
60*9ae013d6SJames Wright 
61*9ae013d6SJames Wright   PetscCall(TSGetDM(ts, &dm));
62*9ae013d6SJames Wright   PetscCall(DMGetOutputDM(dm, &output_dm));
63*9ae013d6SJames Wright   PetscCall(DMGetGlobalVector(output_dm, &Q_output));
64*9ae013d6SJames Wright 
65*9ae013d6SJames Wright   PetscCall(UpdateBoundaryValues(honee, honee->Q_loc, solution_time));
66*9ae013d6SJames Wright   PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, honee->Q_loc));
67*9ae013d6SJames Wright   PetscCall(DMLocalToGlobal(output_dm, honee->Q_loc, INSERT_VALUES, Q_output));
68*9ae013d6SJames Wright 
69*9ae013d6SJames Wright   {  // -- Send solution data to SmartSim
70*9ae013d6SJames Wright     char   array_key[PETSC_MAX_PATH_LEN];
71*9ae013d6SJames Wright     size_t array_key_len;
72*9ae013d6SJames Wright 
73*9ae013d6SJames Wright     if (smartsimsol->overwrite_training_data) {
74*9ae013d6SJames Wright       PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.flow_solution", smartsim->rank_id_name));
75*9ae013d6SJames Wright     } else {
76*9ae013d6SJames Wright       PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.flow_solution.%" PetscInt_FMT, smartsim->rank_id_name, step_num));
77*9ae013d6SJames Wright     }
78*9ae013d6SJames Wright     PetscCall(PetscStrlen(array_key, &array_key_len));
79*9ae013d6SJames Wright 
80*9ae013d6SJames Wright     {
81*9ae013d6SJames Wright       const PetscScalar *Q_output_array;
82*9ae013d6SJames Wright       PetscCall(VecGetArrayRead(Q_output, &Q_output_array));
83*9ae013d6SJames Wright       PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Write, 0, 0, 0, 0));
84*9ae013d6SJames Wright       PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)Q_output_array, smartsimsol->local_array_dims, 2,
85*9ae013d6SJames Wright                                      SRTensorTypeDouble, SRMemLayoutContiguous));
86*9ae013d6SJames Wright       PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Write, 0, 0, 0, 0));
87*9ae013d6SJames Wright       PetscCall(VecRestoreArrayRead(Q_output, &Q_output_array));
88*9ae013d6SJames Wright     }
89*9ae013d6SJames Wright   }
90*9ae013d6SJames Wright 
91*9ae013d6SJames Wright   if (rank % smartsim->collocated_database_num_ranks == 0) {
92*9ae013d6SJames Wright     const char tensor_name[] = "flow_solution_step";
93*9ae013d6SJames Wright     size_t     dim_2[1]      = {2};
94*9ae013d6SJames Wright     PetscInt64 step_array[2] = {step_num, step_num};
95*9ae013d6SJames Wright 
96*9ae013d6SJames Wright     PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
97*9ae013d6SJames Wright     PetscCallSmartRedis(
98*9ae013d6SJames Wright         put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
99*9ae013d6SJames Wright     PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
100*9ae013d6SJames Wright   }
101*9ae013d6SJames Wright 
102*9ae013d6SJames Wright   PetscCall(DMRestoreGlobalVector(output_dm, &Q_output));
103*9ae013d6SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
104*9ae013d6SJames Wright }
105