1 // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 // Based on the instructions from https://www.craylabs.org/docs/sr_integration.html and PHASTA implementation 8 9 #include "../../include/smartsim.h" 10 11 #include "../../navierstokes.h" 12 13 PetscErrorCode SmartRedisVerifyPutTensor(void *c_client, const char *name, const size_t name_length) { 14 bool does_exist = true; 15 16 PetscFunctionBeginUser; 17 PetscSmartRedisCall(tensor_exists(c_client, name, name_length, &does_exist)); 18 PetscCheck(does_exist, PETSC_COMM_SELF, -1, "Tensor of name '%s' was not written to the database successfully", name); 19 PetscFunctionReturn(PETSC_SUCCESS); 20 } 21 22 PetscErrorCode SmartSimTrainingSetup(User user) { 23 SmartSimData smartsim = user->smartsim; 24 PetscMPIInt rank; 25 PetscReal checkrun[2] = {1}; 26 size_t dim_2[1] = {2}; 27 PetscInt num_ranks; 28 29 PetscFunctionBeginUser; 30 PetscCallMPI(MPI_Comm_rank(user->comm, &rank)); 31 PetscCallMPI(MPI_Comm_size(user->comm, &num_ranks)); 32 33 if (rank % smartsim->collocated_database_num_ranks == 0) { 34 // -- Send array that communicates when ML is done training 35 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 36 PetscSmartRedisCall(put_tensor(smartsim->client, "check-run", 9, checkrun, dim_2, 1, SRTensorTypeDouble, SRMemLayoutContiguous)); 37 PetscCall(SmartRedisVerifyPutTensor(smartsim->client, "check-run", 9)); 38 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 39 } 40 PetscFunctionReturn(PETSC_SUCCESS); 41 } 42 43 PetscErrorCode SmartSimSetup(User user) { 44 PetscMPIInt rank; 45 size_t rank_id_name_len; 46 PetscInt num_orchestrator_nodes = 1; 47 48 PetscFunctionBeginUser; 49 PetscCall(PetscNew(&user->smartsim)); 50 SmartSimData smartsim = user->smartsim; 51 52 smartsim->collocated_database_num_ranks = 1; 53 PetscOptionsBegin(user->comm, NULL, "Options for SmartSim integration", NULL); 54 PetscCall(PetscOptionsInt("-smartsim_collocated_database_num_ranks", "Number of ranks per collocated database instance", NULL, 55 smartsim->collocated_database_num_ranks, &smartsim->collocated_database_num_ranks, NULL)); 56 PetscOptionsEnd(); 57 58 PetscCall(PetscStrlen(smartsim->rank_id_name, &rank_id_name_len)); 59 // Create prefix to be put on tensor names 60 PetscCallMPI(MPI_Comm_rank(user->comm, &rank)); 61 PetscCall(PetscSNPrintf(smartsim->rank_id_name, sizeof smartsim->rank_id_name, "y.%d", rank)); 62 63 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Init, 0, 0, 0, 0)); 64 PetscSmartRedisCall(SmartRedisCClient(num_orchestrator_nodes != 1, smartsim->rank_id_name, rank_id_name_len, &smartsim->client)); 65 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Init, 0, 0, 0, 0)); 66 67 PetscCall(SmartSimTrainingSetup(user)); 68 PetscFunctionReturn(PETSC_SUCCESS); 69 } 70 71 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim) { 72 PetscFunctionBeginUser; 73 if (!smartsim) PetscFunctionReturn(PETSC_SUCCESS); 74 75 PetscSmartRedisCall(DeleteCClient(&smartsim->client)); 76 PetscCall(PetscFree(smartsim)); 77 PetscFunctionReturn(PETSC_SUCCESS); 78 } 79