xref: /honee/problems/sgs_dd_model.c (revision 9ab09d51ac9d937f9ba51ebe6202d581c3c26d67)
162b7942eSJames Wright // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors.
262b7942eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
362b7942eSJames Wright //
462b7942eSJames Wright // SPDX-License-Identifier: BSD-2-Clause
562b7942eSJames Wright //
662b7942eSJames Wright // This file is part of CEED:  http://github.com/ceed
762b7942eSJames Wright 
862b7942eSJames Wright #include "../qfunctions/sgs_dd_model.h"
962b7942eSJames Wright 
1062b7942eSJames Wright #include <petscdmplex.h>
1162b7942eSJames Wright 
1262b7942eSJames Wright #include "../navierstokes.h"
1362b7942eSJames Wright 
1462b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN
1562b7942eSJames Wright PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
1662b7942eSJames Wright   PetscFunctionBeginUser;
1762b7942eSJames Wright   for (PetscInt i = 0; i < N; i++) {
1862b7942eSJames Wright     for (PetscInt j = 0; j < M; j++) {
1962b7942eSJames Wright       B[j * N + i] = A[i * M + j];
2062b7942eSJames Wright     }
2162b7942eSJames Wright   }
2262b7942eSJames Wright   PetscFunctionReturn(0);
2362b7942eSJames Wright }
2462b7942eSJames Wright 
2562b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct
2662b7942eSJames Wright PetscErrorCode SGS_DD_ModelContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SGS_DDModelContext *psgsdd_ctx) {
2762b7942eSJames Wright   SGS_DDModelContext sgsdd_ctx;
2862b7942eSJames Wright   PetscInt           num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
2962b7942eSJames Wright   char               file_path[PETSC_MAX_PATH_LEN];
3062b7942eSJames Wright   PetscScalar       *temp;
3162b7942eSJames Wright 
3262b7942eSJames Wright   PetscFunctionBeginUser;
3362b7942eSJames Wright   {
3462b7942eSJames Wright     SGS_DDModelContext sgsdd_temp;
3562b7942eSJames Wright     PetscCall(PetscNew(&sgsdd_temp));
3662b7942eSJames Wright     *sgsdd_temp                     = **psgsdd_ctx;
3762b7942eSJames Wright     sgsdd_temp->offsets.bias1       = 0;
3862b7942eSJames Wright     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
3962b7942eSJames Wright     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
4062b7942eSJames Wright     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
4162b7942eSJames Wright     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
4262b7942eSJames Wright     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
4362b7942eSJames Wright     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
4462b7942eSJames Wright     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
4562b7942eSJames Wright     *sgsdd_ctx = *sgsdd_temp;
4662b7942eSJames Wright     PetscCall(PetscFree(sgsdd_temp));
4762b7942eSJames Wright   }
4862b7942eSJames Wright 
4962b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
5062b7942eSJames Wright   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
5162b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
5262b7942eSJames Wright   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
5362b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
5462b7942eSJames Wright   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
5562b7942eSJames Wright 
5662b7942eSJames Wright   {
5762b7942eSJames Wright     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
5862b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
5962b7942eSJames Wright     PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, temp));
6062b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
6162b7942eSJames Wright     PetscCall(PetscFree(temp));
6262b7942eSJames Wright   }
6362b7942eSJames Wright   {
6462b7942eSJames Wright     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
6562b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
6662b7942eSJames Wright     PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, temp));
6762b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
6862b7942eSJames Wright     PetscCall(PetscFree(temp));
6962b7942eSJames Wright   }
7062b7942eSJames Wright 
7162b7942eSJames Wright   PetscCall(PetscFree(*psgsdd_ctx));
7262b7942eSJames Wright   *psgsdd_ctx = sgsdd_ctx;
7362b7942eSJames Wright   PetscFunctionReturn(0);
7462b7942eSJames Wright }
7562b7942eSJames Wright 
7662b7942eSJames Wright PetscErrorCode SGS_DD_ModelSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
7762b7942eSJames Wright   PetscReal          alpha;
7862b7942eSJames Wright   SGS_DDModelContext sgsdd_ctx;
7962b7942eSJames Wright   MPI_Comm           comm                           = user->comm;
8062b7942eSJames Wright   char               sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_data";
8162b7942eSJames Wright   PetscFunctionBeginUser;
8262b7942eSJames Wright 
83*9ab09d51SJames Wright   PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem));
84*9ab09d51SJames Wright 
8562b7942eSJames Wright   PetscCall(PetscNew(&sgsdd_ctx));
8662b7942eSJames Wright 
8762b7942eSJames Wright   PetscOptionsBegin(comm, NULL, "SGS Data-Drive Model Options", NULL);
8862b7942eSJames Wright   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
8962b7942eSJames Wright   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
9062b7942eSJames Wright                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
9162b7942eSJames Wright   PetscOptionsEnd();
9262b7942eSJames Wright 
9362b7942eSJames Wright   sgsdd_ctx->num_layers  = 2;
9462b7942eSJames Wright   sgsdd_ctx->num_inputs  = 6;
9562b7942eSJames Wright   sgsdd_ctx->num_outputs = 6;
9662b7942eSJames Wright   sgsdd_ctx->num_neurons = 20;
9762b7942eSJames Wright   sgsdd_ctx->alpha       = alpha;
9862b7942eSJames Wright 
9962b7942eSJames Wright   // PetscCall(SGS_DD_ModelContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
10062b7942eSJames Wright 
10162b7942eSJames Wright   PetscFunctionReturn(0);
10262b7942eSJames Wright }
103