Skip to content
Merged
Show file tree
Hide file tree
Changes from 5 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
145 changes: 130 additions & 15 deletions src/Finch_Inputs.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -24,6 +24,7 @@
#include <iostream>
#include <string>
#include <unistd.h>
#include <vector>

#include <nlohmann/json.hpp>

Expand Down Expand Up @@ -168,14 +169,18 @@ class Inputs
MPI_Comm_rank( comm, &comm_rank );
MPI_Comm_size( comm, &comm_size );
std::string filename = getFilename( argc, argv );
parseInputFile( comm, filename );
parseInputFile( filename );
calcAuxiliaryProperties( comm );
}
// constructor for coupled run of Finch with ExaCA
Inputs( MPI_Comm comm, std::string filename )
// constructor for coupled run of Finch with ExaCA - inputs potentially
// spread across both files
Inputs( MPI_Comm comm, const std::string filename,
const int input_file_number = 0 )
{
MPI_Comm_rank( comm, &comm_rank );
MPI_Comm_size( comm, &comm_size );
parseInputFile( comm, filename );
parseInputFile( filename, input_file_number );
calcAuxiliaryProperties( comm );
}

void write()
Expand Down Expand Up @@ -264,12 +269,81 @@ class Inputs
return filename_s;
}

void parseInputFile( MPI_Comm comm, const std::string filename )
void parseInputFile( const std::string filename,
const int input_file_number = 0 )
{
readInput( filename );

// Input file is either a Finch input file or an ExaCA input file with a
// Finch object
Info << "Parsing input file " << input_file_number << std::endl;
std::ifstream input_data_stream( filename );
nlohmann::json input_data_raw =
nlohmann::json::parse( input_data_stream );
if ( !input_data_raw.contains( "Finch" ) )
{
// This is a Finch input file and should have all 5 sections
std::vector<bool> found_sections = readSections( input_data_raw );
if ( std::find( found_sections.begin(), found_sections.end(),
false ) != found_sections.end() )
throw std::runtime_error(
"Error: Missing top-level sections of Finch input file, "
"see README for proper input file format" );
}
else
{
// This is an ExaCA input file - get data from the Finch object
nlohmann::json top_level_input_data = input_data_raw["Finch"];
// Sections to parse
std::vector<std::string> input_file_sections = {
"time", "space", "properties", "source", "sampling" };
const int num_inp_file_sections = input_file_sections.size();
if ( !top_level_input_data.contains( "layers" ) )
{
// All inputs should be present at the top level of the file
std::vector<bool> found_sections =
readSections( top_level_input_data );
if ( std::find( found_sections.begin(), found_sections.end(),
false ) != found_sections.end() )
throw std::runtime_error(
"Error: Missing top-level sections of Finch input "
"file, see README for proper input file format" );
}
else
{
// Check top level for inputs and parse
std::vector<bool> found_sections_top_level =
readSections( top_level_input_data );
// Check individual layer level for inputs and parse
std::vector<bool> found_sections_layer_level = readSections(
top_level_input_data["layers"][input_file_number] );
// Ensure each input was present at least once, warn about
// redundant inputs (layer level takes priority)
for ( int n = 0; n < num_inp_file_sections; n++ )
{
if ( ( found_sections_top_level[n] ) &&
( found_sections_layer_level[n] ) )
{
Info << "Warning: Finch input object "
<< input_file_sections[n]
<< " has multiple values given; values from "
"`layers` object will be used"
<< std::endl;
}
else if ( ( !found_sections_top_level[n] ) &&
( !found_sections_layer_level[n] ) )
{
std::string err_message = "Error: Finch input object " +
input_file_sections[n] +
" was not found";
throw std::runtime_error( err_message );
}
}
}
}
write();
}

void calcAuxiliaryProperties( MPI_Comm comm )
{
// create auxiliary properties
properties.thermal_diffusivity =
( properties.thermal_conductivity ) /
Expand All @@ -292,30 +366,62 @@ class Inputs
time_monitor = TimeMonitor( comm, time );
}

void readInput( const std::string filename )
// Calls other read input functions to initialize variables, returning a
// list of which sections were found
std::vector<bool> readSections( nlohmann::json db )
{
// parse input file
std::ifstream db_stream( filename );
nlohmann::json db = nlohmann::json::parse( db_stream );
std::vector<bool> found_sections( 5, false );
Comment thread
streeve marked this conversation as resolved.
if ( db.contains( "time" ) )
{
readInput_Time( db );
found_sections[0] = true;
Comment thread
streeve marked this conversation as resolved.
}
if ( db.contains( "space" ) )
{
readInput_Space( db );
found_sections[1] = true;
}
if ( db.contains( "properties" ) )
{
readInput_Properties( db );
found_sections[2] = true;
}
if ( db.contains( "source" ) )
{
readInput_Source( db );
found_sections[3] = true;
}
if ( db.contains( "sampling" ) )
{
readInput_Sampling( db );
found_sections[4] = true;
}
return found_sections;
}

void readInput_Time( nlohmann::json db )
Comment thread
streeve marked this conversation as resolved.
Outdated
{
// Read time components
time.Co = db["time"]["Co"];
time.start_time = db["time"]["start_time"];
time.end_time = db["time"]["end_time"];
time.output.total_steps = db["time"]["total_output_steps"];
time.monitor.total_steps = db["time"]["total_monitor_steps"];
}

void readInput_Space( nlohmann::json db )
{
// Read space components
space.initial_temperature = db["space"]["initial_temperature"];
space.cell_size = db["space"]["cell_size"];
space.global_low_corner = db["space"]["global_low_corner"];
space.global_high_corner = db["space"]["global_high_corner"];

/*
Default block partitioner. This relies on MPI_Cart_create to
balance the number of ranks in each direction. This partitioning
is best only in the global mesh is a uniform cube.
*/
Default block partitioner. This relies on MPI_Cart_create to
balance the number of ranks in each direction. This partitioning
is best only in the global mesh is a uniform cube.
*/
std::array<int, 3> default_ranks_per_dim = { 0, 0, 0 };

std::array<int, 3> ranks_per_dim = default_ranks_per_dim;
Expand All @@ -331,7 +437,10 @@ class Inputs
}

space.ranks_per_dim = ranks_per_dim;
}

void readInput_Properties( nlohmann::json db )
{
// Read properties components
properties.density = db["properties"]["density"];
properties.specific_heat = db["properties"]["specific_heat"];
Expand All @@ -340,7 +449,10 @@ class Inputs
properties.latent_heat = db["properties"]["latent_heat"];
properties.solidus = db["properties"]["solidus"];
properties.liquidus = db["properties"]["liquidus"];
}

void readInput_Source( nlohmann::json db )
{
// Read heat source components
source.absorption = db["source"]["absorption"];
source.two_sigma = db["source"]["two_sigma"];
Expand All @@ -350,7 +462,10 @@ class Inputs
source.two_sigma[2] = fabs( source.two_sigma[2] );

source.scan_path_file = db["source"]["scan_path_file"];
}

void readInput_Sampling( nlohmann::json db )
{
// Read sampling components
sampling.enabled = false;
if ( db.contains( "sampling" ) )
Expand Down
44 changes: 44 additions & 0 deletions src/Finch_Run.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -107,6 +107,50 @@ class Layer
}

auto getSolidificationData() { return solidification_data_.get(); }
Comment thread
streeve marked this conversation as resolved.
// Append next layer's solidification data to input_solidification_data
void appendSolidificationData(
Kokkos::View<double**, Kokkos::LayoutLeft, Kokkos::HostSpace>&
input_solidification_data,
std::vector<int>& first_value_finch, std::vector<int>& last_value_finch,
Comment thread
streeve marked this conversation as resolved.
int finch_file_num, const int num_finch_simulations )
{
// Time-temperature history from the Finch simulation performed for this
// layer
auto new_layer_data = solidification_data_.get();
// Number of events and components in new layer time-temperature history
const int events_this_layer = new_layer_data.extent( 0 );
const int n_cmpts = new_layer_data.extent( 1 );
// Number of events in currently stored time-temperature history-
// first_value_finch provides offset from data stored for previous
// layers if performing more than 1 finch simulation at a time
int events_prev_layers;
if ( ( finch_file_num == 0 ) || ( num_finch_simulations == 1 ) )
{
first_value_finch[finch_file_num] = 0;
events_prev_layers = 0;
}
else
{
first_value_finch[finch_file_num] =
last_value_finch[finch_file_num - 1];
events_prev_layers = input_solidification_data.extent( 0 );
}
// Resize input_solidification_data to accommodate both any events
// stored from previous layers and the events calculated from simulation
// of this layer
Kokkos::resize( input_solidification_data,
events_prev_layers + events_this_layer, n_cmpts );
// Copy this layer's data into the return view
for ( int i = 0; i < events_this_layer; i++ )
for ( int j = 0; j < n_cmpts; j++ )
input_solidification_data(
first_value_finch[finch_file_num] + i, j ) =
new_layer_data( i, j );
// Set last_value_finch to bound the indices with time-temperature
// history data for this layer
last_value_finch[finch_file_num] =
events_prev_layers + events_this_layer;
}

auto writeSolidificationData( Sampling sampling_inputs, MPI_Comm comm )
{
Expand Down
Loading