19ae013d6SJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2025, HONEE contributors. 29ae013d6SJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 39ae013d6SJames Wright 49ae013d6SJames Wright #include <navierstokes.h> 59ae013d6SJames Wright #include <smartsim-impl.h> 69ae013d6SJames Wright 79ae013d6SJames Wright typedef struct { 89ae013d6SJames Wright PetscInt write_data_interval; 99ae013d6SJames Wright size_t local_array_dims[2]; 109ae013d6SJames Wright PetscBool overwrite_training_data; 119ae013d6SJames Wright } *SmartSimSolutionData; 129ae013d6SJames Wright 139ae013d6SJames Wright PetscErrorCode TSMonitor_SmartSimSolutionSetup(TS ts, PetscViewerAndFormat *ctx) { 149ae013d6SJames Wright SmartSimSolutionData smartsimsol; 159ae013d6SJames Wright Honee honee; 169ae013d6SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 179ae013d6SJames Wright 189ae013d6SJames Wright PetscFunctionBeginUser; 199ae013d6SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 209ae013d6SJames Wright PetscCall(PetscNew(&smartsimsol)); 219ae013d6SJames Wright smartsimsol->overwrite_training_data = PETSC_TRUE; 229ae013d6SJames Wright 239ae013d6SJames Wright PetscOptionsBegin(comm, NULL, "SmartSim Solution Writing", NULL); 249ae013d6SJames Wright PetscCall(PetscOptionsBool("-ts_monitor_smartsim_solution_overwrite_data", "Overwrite old solution data in the database", NULL, 259ae013d6SJames Wright smartsimsol->overwrite_training_data, &smartsimsol->overwrite_training_data, NULL)); 269ae013d6SJames Wright PetscOptionsEnd(); 279ae013d6SJames Wright 289ae013d6SJames Wright { // Get solution vector size 299ae013d6SJames Wright PetscSection output_section; 309ae013d6SJames Wright PetscInt num_dofs, num_comps; 319ae013d6SJames Wright DM dm, output_dm; 329ae013d6SJames Wright 339ae013d6SJames Wright PetscCall(TSGetDM(ts, &dm)); 349ae013d6SJames Wright PetscCall(DMGetOutputDM(dm, &output_dm)); 359ae013d6SJames Wright 369ae013d6SJames Wright PetscCall(DMGetGlobalVectorInfo(output_dm, &num_dofs, NULL, NULL)); 379ae013d6SJames Wright PetscCall(DMGetGlobalSection(output_dm, &output_section)); 389ae013d6SJames Wright PetscCall(PetscSectionGetFieldComponents(output_section, 0, &num_comps)); 399ae013d6SJames Wright smartsimsol->local_array_dims[0] = num_dofs / num_comps; 409ae013d6SJames Wright smartsimsol->local_array_dims[1] = num_comps; 419ae013d6SJames Wright } 429ae013d6SJames Wright ctx->data = smartsimsol; 439ae013d6SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 449ae013d6SJames Wright } 459ae013d6SJames Wright 469ae013d6SJames Wright PetscErrorCode TSMonitor_SmartSimSolution(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) { 479ae013d6SJames Wright Honee honee; 489ae013d6SJames Wright SmartSimData smartsim; 499ae013d6SJames Wright Vec Q_output; 509ae013d6SJames Wright PetscMPIInt rank; 519ae013d6SJames Wright DM dm, output_dm; 529ae013d6SJames Wright 539ae013d6SJames Wright PetscFunctionBeginUser; 549ae013d6SJames Wright SmartSimSolutionData smartsimsol = ctx->data; 559ae013d6SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 569ae013d6SJames Wright PetscCall(HoneeGetSmartSimData(honee, &smartsim)); 579ae013d6SJames Wright PetscCallMPI(MPI_Comm_rank(honee->comm, &rank)); 589ae013d6SJames Wright 599ae013d6SJames Wright if (step_num % ctx->view_interval != 0) PetscFunctionReturn(PETSC_SUCCESS); 609ae013d6SJames Wright 619ae013d6SJames Wright PetscCall(TSGetDM(ts, &dm)); 629ae013d6SJames Wright PetscCall(DMGetOutputDM(dm, &output_dm)); 639ae013d6SJames Wright PetscCall(DMGetGlobalVector(output_dm, &Q_output)); 649ae013d6SJames Wright 659ae013d6SJames Wright PetscCall(UpdateBoundaryValues(honee, honee->Q_loc, solution_time)); 669ae013d6SJames Wright PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, honee->Q_loc)); 679ae013d6SJames Wright PetscCall(DMLocalToGlobal(output_dm, honee->Q_loc, INSERT_VALUES, Q_output)); 689ae013d6SJames Wright 699ae013d6SJames Wright { // -- Send solution data to SmartSim 709ae013d6SJames Wright char array_key[PETSC_MAX_PATH_LEN]; 719ae013d6SJames Wright size_t array_key_len; 72*1fc84d60SJames Wright void *dataset; 739ae013d6SJames Wright 749ae013d6SJames Wright if (smartsimsol->overwrite_training_data) { 759ae013d6SJames Wright PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.flow_solution", smartsim->rank_id_name)); 769ae013d6SJames Wright } else { 779ae013d6SJames Wright PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.flow_solution.%" PetscInt_FMT, smartsim->rank_id_name, step_num)); 789ae013d6SJames Wright } 799ae013d6SJames Wright PetscCall(PetscStrlen(array_key, &array_key_len)); 80*1fc84d60SJames Wright PetscCallSmartRedis(CDataSet(array_key, array_key_len, &dataset)); 819ae013d6SJames Wright 829ae013d6SJames Wright { 839ae013d6SJames Wright const PetscScalar *Q_output_array; 849ae013d6SJames Wright PetscCall(VecGetArrayRead(Q_output, &Q_output_array)); 859ae013d6SJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Write, 0, 0, 0, 0)); 86*1fc84d60SJames Wright PetscCallSmartRedis( 87*1fc84d60SJames Wright add_tensor(dataset, "solution", 8, (void *)Q_output_array, smartsimsol->local_array_dims, 2, SRTensorTypeDouble, SRMemLayoutContiguous)); 889ae013d6SJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Write, 0, 0, 0, 0)); 899ae013d6SJames Wright PetscCall(VecRestoreArrayRead(Q_output, &Q_output_array)); 909ae013d6SJames Wright } 919ae013d6SJames Wright 92*1fc84d60SJames Wright PetscCallSmartRedis(add_meta_scalar(dataset, "step", 4, (void *)&step_num, SRMetadataTypePetscInt)); 93*1fc84d60SJames Wright PetscCallSmartRedis(add_meta_scalar(dataset, "time", 4, (void *)&solution_time, SRMetadataTypeDouble)); 94*1fc84d60SJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Write, 0, 0, 0, 0)); 95*1fc84d60SJames Wright PetscCallSmartRedis(put_dataset(smartsim->client, dataset)); 96*1fc84d60SJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Write, 0, 0, 0, 0)); 97*1fc84d60SJames Wright PetscCallSmartRedis(DeallocateDataSet(&dataset)); 989ae013d6SJames Wright } 999ae013d6SJames Wright 1009ae013d6SJames Wright PetscCall(DMRestoreGlobalVector(output_dm, &Q_output)); 1019ae013d6SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1029ae013d6SJames Wright } 103