Skip to content

Commit

Permalink
Merge pull request #551 from PrincetonUniversity/issue-550
Browse files Browse the repository at this point in the history
Merge kernel and property writer
  • Loading branch information
icui authored Mar 5, 2025
2 parents ac84de2 + ae21b6f commit 7923695
Show file tree
Hide file tree
Showing 6 changed files with 159 additions and 217 deletions.
14 changes: 14 additions & 0 deletions include/IO/impl/medium_writer.hpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,14 @@
#pragma once

namespace specfem {
namespace IO {
namespace impl {
template <typename OutputLibrary, typename ContainerType>
void write_container(const std::string &output_folder,
const std::string &output_namespace,
const specfem::compute::mesh &mesh,
const specfem::compute::element_types &element_types,
ContainerType &container);
} // namespace impl
} // namespace IO
} // namespace specfem
118 changes: 118 additions & 0 deletions include/IO/impl/medium_writer.tpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,118 @@
#pragma once

#include "IO/impl/medium_writer.hpp"
#include "compute/assembly/assembly.hpp"
#include "enumerations/dimension.hpp"
#include "enumerations/medium.hpp"
#include "kokkos_abstractions.h"
#include <Kokkos_Core.hpp>

template <typename OutputLibrary, typename ContainerType>
void specfem::IO::impl::write_container(
const std::string &output_folder, const std::string &output_namespace,
const specfem::compute::mesh &mesh,
const specfem::compute::element_types &element_types,
ContainerType &container) {
using DomainView =
Kokkos::View<type_real ***, Kokkos::LayoutLeft, Kokkos::HostSpace>;

container.copy_to_host();

typename OutputLibrary::File file(output_folder + "/" + output_namespace);

const int nspec = mesh.points.nspec;
const int ngllz = mesh.points.ngllz;
const int ngllx = mesh.points.ngllx;

int n_elastic_isotropic;
int n_elastic_anisotropic;
int n_acoustic;

{
typename OutputLibrary::Group elastic =
file.createGroup("/ElasticIsotropic");

const auto element_indices = element_types.get_elements_on_host(
specfem::element::medium_tag::elastic,
specfem::element::property_tag::isotropic);
n_elastic_isotropic = element_indices.size();

DomainView x("xcoordinates_elastic_isotropic", n_elastic_isotropic, ngllz,
ngllx);
DomainView z("zcoordinates_elastic_isotropic", n_elastic_isotropic, ngllz,
ngllx);

for (int i = 0; i < n_elastic_isotropic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

elastic.createDataset("X", x).write();
elastic.createDataset("Z", z).write();
elastic.createDataset("data", container.elastic_isotropic.h_data).write();
}

{
typename OutputLibrary::Group elastic =
file.createGroup("/ElasticAnisotropic");

const auto element_indices = element_types.get_elements_on_host(
specfem::element::medium_tag::elastic,
specfem::element::property_tag::anisotropic);
n_elastic_anisotropic = element_indices.size();

DomainView x("xcoordinates_elastic_anisotropic", n_elastic_anisotropic,
ngllz, ngllx);
DomainView z("zcoordinates_elastic_anisotropic", n_elastic_anisotropic,
ngllz, ngllx);

for (int i = 0; i < n_elastic_anisotropic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

elastic.createDataset("X", x).write();
elastic.createDataset("Z", z).write();
elastic.createDataset("data", container.elastic_anisotropic.h_data).write();
}

{
typename OutputLibrary::Group acoustic = file.createGroup("/Acoustic");

const auto element_indices = element_types.get_elements_on_host(
specfem::element::medium_tag::acoustic);
n_acoustic = element_indices.size();

DomainView x("xcoordinates_acoustic", n_acoustic, ngllz, ngllx);
DomainView z("zcoordinates_acoustic", n_acoustic, ngllz, ngllx);

for (int i = 0; i < n_acoustic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

acoustic.createDataset("X", x).write();
acoustic.createDataset("Z", z).write();
acoustic.createDataset("data", container.acoustic_isotropic.h_data).write();
}

assert(n_elastic_isotropic + n_elastic_anisotropic + n_acoustic == nspec);

std::cout << output_namespace << " written to " << output_folder << "/"
<< output_namespace << std::endl;
}
119 changes: 11 additions & 108 deletions include/IO/kernel/writer.tpp
Original file line number Diff line number Diff line change
@@ -1,116 +1,19 @@
#pragma once

#include "compute/assembly/assembly.hpp"
#include "enumerations/dimension.hpp"
#include "enumerations/medium.hpp"
#include "kokkos_abstractions.h"
#include "point/kernels.hpp"
#include "IO/kernel/writer.hpp"
#include <Kokkos_Core.hpp>
#include "IO/property/writer.hpp"
#include "IO/impl/medium_writer.hpp"

namespace specfem {
namespace IO {

template <typename OutputLibrary>
specfem::IO::kernel_writer<OutputLibrary>::kernel_writer(const std::string output_folder)
kernel_writer<OutputLibrary>::kernel_writer(const std::string output_folder)
: output_folder(output_folder) {}

template <typename OutputLibrary>
void specfem::IO::kernel_writer<OutputLibrary>::write(specfem::compute::assembly &assembly) {
const auto &mesh = assembly.mesh;
auto &element_types = assembly.element_types;
auto &kernels = assembly.kernels;

using DomainView =
Kokkos::View<type_real ***, Kokkos::LayoutLeft, Kokkos::HostSpace>;

kernels.copy_to_host();

typename OutputLibrary::File file(output_folder + "/Kernels");

const int nspec = mesh.points.nspec;
const int ngllz = mesh.points.ngllz;
const int ngllx = mesh.points.ngllx;

int n_elastic_isotropic;
int n_elastic_anisotropic;
int n_acoustic;

{
typename OutputLibrary::Group elastic = file.createGroup("/ElasticIsotropic");

const auto element_indices = element_types.get_elements_on_host(
specfem::element::medium_tag::elastic,
specfem::element::property_tag::isotropic);
n_elastic_isotropic = element_indices.size();

DomainView x("xcoordinates_elastic_isotropic", n_elastic_isotropic, ngllz, ngllx);
DomainView z("zcoordinates_elastic_isotropic", n_elastic_isotropic, ngllz, ngllx);

for (int i = 0; i < n_elastic_isotropic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

elastic.createDataset("X", x).write();
elastic.createDataset("Z", z).write();
elastic.createDataset("data", kernels.elastic_isotropic.h_data).write();
}

{
typename OutputLibrary::Group elastic = file.createGroup("/ElasticAnisotropic");

const auto element_indices = element_types.get_elements_on_host(
specfem::element::medium_tag::elastic,
specfem::element::property_tag::anisotropic);
n_elastic_anisotropic = element_indices.size();

DomainView x("xcoordinates_elastic_anisotropic", n_elastic_anisotropic, ngllz, ngllx);
DomainView z("zcoordinates_elastic_anisotropic", n_elastic_anisotropic, ngllz, ngllx);

for (int i = 0; i < n_elastic_anisotropic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

elastic.createDataset("X", x).write();
elastic.createDataset("Z", z).write();
elastic.createDataset("data", kernels.elastic_anisotropic.h_data).write();
}

{
typename OutputLibrary::Group acoustic = file.createGroup("/Acoustic");

const auto element_indices = element_types.get_elements_on_host(specfem::element::medium_tag::acoustic);
n_acoustic = element_indices.size();

DomainView x("xcoordinates_acoustic", n_acoustic, ngllz, ngllx);
DomainView z("zcoordinates_acoustic", n_acoustic, ngllz, ngllx);

for (int i = 0; i < n_acoustic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

acoustic.createDataset("X", x).write();
acoustic.createDataset("Z", z).write();
acoustic.createDataset("data", kernels.acoustic_isotropic.h_data).write();
}

assert(n_elastic_isotropic + n_elastic_anisotropic + n_acoustic == nspec);

std::cout << "Kernels written to " << output_folder << "/Kernels"
<< std::endl;
void kernel_writer<OutputLibrary>::write(specfem::compute::assembly &assembly) {
impl::write_container<OutputLibrary>(output_folder, "Kernels", assembly.mesh, assembly.element_types, assembly.kernels);
}

} // namespace IO
} // namespace specfem
117 changes: 10 additions & 107 deletions include/IO/property/writer.tpp
Original file line number Diff line number Diff line change
@@ -1,116 +1,19 @@
#pragma once

#include "compute/assembly/assembly.hpp"
#include "enumerations/dimension.hpp"
#include "enumerations/medium.hpp"
#include "kokkos_abstractions.h"
#include "point/properties.hpp"
#include "IO/property/writer.hpp"
#include <Kokkos_Core.hpp>
#include "IO/impl/medium_writer.hpp"

namespace specfem {
namespace IO {

template <typename OutputLibrary>
specfem::IO::property_writer<OutputLibrary>::property_writer(const std::string output_folder)
property_writer<OutputLibrary>::property_writer(const std::string output_folder)
: output_folder(output_folder) {}

template <typename OutputLibrary>
void specfem::IO::property_writer<OutputLibrary>::write(specfem::compute::assembly &assembly) {
const auto &mesh = assembly.mesh;
auto &element_types = assembly.element_types;
auto &properties = assembly.properties;

using DomainView =
Kokkos::View<type_real ***, Kokkos::LayoutLeft, Kokkos::HostSpace>;

properties.copy_to_host();

typename OutputLibrary::File file(output_folder + "/Properties");

const int nspec = mesh.points.nspec;
const int ngllz = mesh.points.ngllz;
const int ngllx = mesh.points.ngllx;

int n_elastic_isotropic;
int n_elastic_anisotropic;
int n_acoustic;

{
typename OutputLibrary::Group elastic = file.createGroup("/ElasticIsotropic");

const auto element_indices = element_types.get_elements_on_host(
specfem::element::medium_tag::elastic,
specfem::element::property_tag::isotropic);
n_elastic_isotropic = element_indices.size();

DomainView x("xcoordinates_elastic_isotropic", n_elastic_isotropic, ngllz, ngllx);
DomainView z("zcoordinates_elastic_isotropic", n_elastic_isotropic, ngllz, ngllx);

for (int i = 0; i < n_elastic_isotropic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

elastic.createDataset("X", x).write();
elastic.createDataset("Z", z).write();
elastic.createDataset("data", properties.elastic_isotropic.h_data).write();
}

{
typename OutputLibrary::Group elastic = file.createGroup("/ElasticAnisotropic");

const auto element_indices = element_types.get_elements_on_host(
specfem::element::medium_tag::elastic,
specfem::element::property_tag::anisotropic);
n_elastic_anisotropic = element_indices.size();

DomainView x("xcoordinates_elastic_anisotropic", n_elastic_anisotropic, ngllz, ngllx);
DomainView z("zcoordinates_elastic_anisotropic", n_elastic_anisotropic, ngllz, ngllx);

for (int i = 0; i < n_elastic_anisotropic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

elastic.createDataset("X", x).write();
elastic.createDataset("Z", z).write();
elastic.createDataset("data", properties.elastic_anisotropic.h_data).write();
}

{
typename OutputLibrary::Group acoustic = file.createGroup("/Acoustic");

const auto element_indices = element_types.get_elements_on_host(specfem::element::medium_tag::acoustic);
n_acoustic = element_indices.size();

DomainView x("xcoordinates_acoustic", n_acoustic, ngllz, ngllx);
DomainView z("zcoordinates_acoustic", n_acoustic, ngllz, ngllx);

for (int i = 0; i < n_acoustic; i++) {
const int ispec = element_indices(i);
for (int iz = 0; iz < ngllz; iz++) {
for (int ix = 0; ix < ngllx; ix++) {
x(i, iz, ix) = mesh.points.h_coord(0, ispec, iz, ix);
z(i, iz, ix) = mesh.points.h_coord(1, ispec, iz, ix);
}
}
}

acoustic.createDataset("X", x).write();
acoustic.createDataset("Z", z).write();
acoustic.createDataset("data", properties.acoustic_isotropic.h_data).write();
}

assert(n_elastic_isotropic + n_elastic_anisotropic + n_acoustic == nspec);

std::cout << "Properties written to " << output_folder << "/Properties"
<< std::endl;
void property_writer<OutputLibrary>::write(specfem::compute::assembly &assembly) {
impl::write_container<OutputLibrary>(output_folder, "Properties", assembly.mesh, assembly.element_types, assembly.properties);
}

} // namespace IO
} // namespace specfem
Loading

0 comments on commit 7923695

Please sign in to comment.