have first cut of just pair-dpd

This commit is contained in:
Matt Bettencourt
2022-03-30 16:22:26 +02:00
parent ddf93eb7ba
commit 5003c35963
5 changed files with 585 additions and 4 deletions

View File

@ -46,6 +46,8 @@ PairDPD::PairDPD(LAMMPS *lmp) : Pair(lmp)
PairDPD::~PairDPD()
{
if (copymode) return;
if (allocated) {
memory->destroy(setflag);
memory->destroy(cutsq);

View File

@ -28,11 +28,11 @@ class PairDPD : public Pair {
public:
PairDPD(class LAMMPS *);
~PairDPD() override;
void compute(int, int) override;
virtual void compute(int, int) override;
void settings(int, char **) override;
void coeff(int, char **) override;
void init_style() override;
double init_one(int, int) override;
virtual void init_style() override;
virtual double init_one(int, int) override;
void write_restart(FILE *) override;
void read_restart(FILE *) override;
void write_restart_settings(FILE *) override;
@ -49,7 +49,7 @@ class PairDPD : public Pair {
double **sigma;
class RanMars *random;
void allocate();
virtual void allocate();
};
} // namespace LAMMPS_NS

View File

@ -231,6 +231,8 @@ action pair_coul_long_kokkos.h pair_coul_long.h
action pair_coul_wolf_kokkos.cpp
action pair_coul_wolf_kokkos.h
action pair_dpd_ext_kokkos.cpp pair_dpd_ext.cpp
action pair_dpd_kokkos.h pair_dpd.h
action pair_dpd_kokkos.cpp pair_dpd.cpp
action pair_dpd_ext_kokkos.h pair_dpd_ext.h
action pair_dpd_fdt_energy_kokkos.cpp pair_dpd_fdt_energy.cpp
action pair_dpd_fdt_energy_kokkos.h pair_dpd_fdt_energy.h

View File

@ -0,0 +1,451 @@
// clang-format off
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
https://www.lammps.org/, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
/* ----------------------------------------------------------------------
Contributing author: Matt Bettencourt (NVIDIA)
------------------------------------------------------------------------- */
#include "pair_dpd_kokkos.h"
#include "atom.h"
#include "atom_kokkos.h"
#include "memory_kokkos.h"
#include "comm.h"
#include "error.h"
#include "force.h"
#include "memory.h"
#include "neigh_list.h"
#include "neigh_request.h"
#include "neighbor.h"
#include "random_mars.h"
#include "update.h"
#include "atom_masks.h"
#include "kokkos.h"
#include <cmath>
using namespace LAMMPS_NS;
#define EPSILON 1.0e-10
template<class DeviceType>
PairDPDKokkos<DeviceType>::PairDPDKokkos(class LAMMPS *lmp) :
PairDPD(lmp) ,
#ifdef DPD_USE_RAN_MARS
rand_pool(0 /* unused */, lmp)
#else
rand_pool()
#endif
{
kokkosable = 1;
atomKK = (AtomKokkos *) atom;
execution_space = ExecutionSpaceFromDevice<DeviceType>::space;
datamask_read = X_MASK | V_MASK | F_MASK | TYPE_MASK | ENERGY_MASK | VIRIAL_MASK;
datamask_modify = F_MASK | ENERGY_MASK | VIRIAL_MASK;
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
PairDPDKokkos<DeviceType>::~PairDPDKokkos() {
if (copymode) return;
#ifdef DPD_USE_RAN_MARS
rand_pool.destroy();
#endif
memoryKK->destroy_kokkos(k_eatom,eatom);
memoryKK->destroy_kokkos(k_vatom,vatom);
memoryKK->destroy_kokkos(k_cutsq,cutsq);
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
void PairDPDKokkos<DeviceType>::init_style()
{
PairDPD::init_style();
#ifdef DPD_USE_RAN_MARS
rand_pool.init(random,seed);
#else
typedef Kokkos::Experimental::UniqueToken<
DeviceType, Kokkos::Experimental::UniqueTokenScope::Global> unique_token_type;
unique_token_type unique_token;
rand_pool.init(seed + comm->me,unique_token.size());
#endif
neighflag = lmp->kokkos->neighflag;
if (force->newton_pair == 0 || neighflag == FULL )
error->all(FLERR,"Must use half neighbor list style and newton on with pair dpd/kk");
auto request = neighbor->find_request(this);
request->set_kokkos_host(std::is_same<DeviceType,LMPHostType>::value &&
!std::is_same<DeviceType,LMPDeviceType>::value);
request->set_kokkos_device(std::is_same<DeviceType,LMPDeviceType>::value);
if (neighflag == FULL)
request->enable_full();
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
void PairDPDKokkos<DeviceType>::compute(int eflagin, int vflagin)
{
eflag = eflagin; vflag = vflagin;
if (neighflag == FULL) no_virial_fdotr_compute = 1;
ev_init(eflag,vflag);
if (eflag_atom) {
memoryKK->destroy_kokkos(k_eatom,eatom);
memoryKK->create_kokkos(k_eatom,eatom,maxeatom,"pair:eatom");
d_eatom = k_eatom.template view<DeviceType>();
}
if (vflag_atom) {
memoryKK->destroy_kokkos(k_vatom,vatom);
memoryKK->create_kokkos(k_vatom,vatom,maxvatom,"pair:vatom");
d_vatom = k_vatom.template view<DeviceType>();
}
atomKK->sync(execution_space,X_MASK | V_MASK | F_MASK | TYPE_MASK | ENERGY_MASK | VIRIAL_MASK);
x = atomKK->k_x.view<DeviceType>();
v = atomKK->k_v.view<DeviceType>();
f = atomKK->k_f.view<DeviceType>();
type = atomKK->k_type.view<DeviceType>();
k_cutsq.template sync<DeviceType>();
k_params.template sync<DeviceType>();
special_lj[0] = force->special_lj[0];
special_lj[1] = force->special_lj[1];
special_lj[2] = force->special_lj[2];
special_lj[3] = force->special_lj[3];
nlocal = atom->nlocal;
newton_pair = force->newton_pair;
dtinvsqrt = 1.0/sqrt(update->dt);
NeighListKokkos<DeviceType>* k_list = static_cast<NeighListKokkos<DeviceType>*>(list);
d_numneigh = k_list->d_numneigh;
d_neighbors = k_list->d_neighbors;
d_ilist = k_list->d_ilist;
// loop over neighbors of my atoms
int inum = list->inum;
EV_FLOAT ev;
copymode = 1;
if (neighflag == HALF) {
if (newton_pair) {
if (evflag) Kokkos::parallel_reduce(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<HALF,1,1> >(0,inum),*this,ev);
else Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<HALF,1,0> >(0,inum),*this);
} else {
if (evflag) Kokkos::parallel_reduce(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<HALF,0,1> >(0,inum),*this,ev);
else Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<HALF,0,0> >(0,inum),*this);
}
} else if (neighflag == HALFTHREAD) {
if (newton_pair) {
if (evflag) Kokkos::parallel_reduce(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<HALFTHREAD,1,1> >(0,inum),*this,ev);
else Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<HALFTHREAD,1,0> >(0,inum),*this);
} else {
if (evflag) Kokkos::parallel_reduce(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<HALFTHREAD,0,1> >(0,inum),*this,ev);
else Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<HALFTHREAD,0,0> >(0,inum),*this);
}
} else if (neighflag == FULL) {
if (newton_pair) {
if (evflag) Kokkos::parallel_reduce(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<FULL,1,1> >(0,inum),*this,ev);
else Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<FULL,1,0> >(0,inum),*this);
} else {
if (evflag) Kokkos::parallel_reduce(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<FULL,0,1> >(0,inum),*this,ev);
else Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType, TagDPDKokkos<FULL,0,0> >(0,inum),*this);
}
}
if (eflag_global) eng_vdwl += ev.evdwl;
if (vflag_global) {
virial[0] += ev.v[0];
virial[1] += ev.v[1];
virial[2] += ev.v[2];
virial[3] += ev.v[3];
virial[4] += ev.v[4];
virial[5] += ev.v[5];
}
if (vflag_fdotr) pair_virial_fdotr_compute(this);
if (eflag_atom) {
k_eatom.template modify<DeviceType>();
k_eatom.template sync<LMPHostType>();
}
if (vflag_atom) {
k_vatom.template modify<DeviceType>();
k_vatom.template sync<LMPHostType>();
}
copymode = 0;
if (evflag) atomKK->modified(execution_space,F_MASK | ENERGY_MASK | VIRIAL_MASK);
else atomKK->modified(execution_space,F_MASK);
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
template<int NEIGHFLAG, int NEWTON_PAIR, int EVFLAG>
KOKKOS_INLINE_FUNCTION
void PairDPDKokkos<DeviceType>::operator() (TagDPDKokkos<NEIGHFLAG,NEWTON_PAIR,EVFLAG>, const int &ii) const {
EV_FLOAT ev;
this->template operator()<NEIGHFLAG,NEWTON_PAIR,EVFLAG>(TagDPDKokkos<NEIGHFLAG,NEWTON_PAIR,EVFLAG>(), ii, ev);
}
template<class DeviceType>
template<int NEIGHFLAG, int NEWTON_PAIR, int EVFLAG>
KOKKOS_INLINE_FUNCTION
void PairDPDKokkos<DeviceType>::operator() (TagDPDKokkos<NEIGHFLAG,NEWTON_PAIR,EVFLAG>, const int &ii, EV_FLOAT &ev) const {
Kokkos::View<F_FLOAT*[3], typename DAT::t_f_array::array_layout,typename KKDevice<DeviceType>::value,Kokkos::MemoryTraits<AtomicF<NEIGHFLAG>::value> > a_f = f;
int i,j,jj,jnum,itype,jtype;
double xtmp,ytmp,ztmp,delx,dely,delz,fpair;
double vxtmp,vytmp,vztmp,delvx,delvy,delvz;
double rsq,r,rinv,dot,wd,randnum,factor_dpd;
double fx = 0,fy = 0,fz = 0;
double evdwl = 0;
i = d_ilist[ii];
xtmp = x(i,0);
ytmp = x(i,1);
ztmp = x(i,2);
vxtmp = v(i,0);
vytmp = v(i,1);
vztmp = v(i,2);
itype = type(i);
jnum = d_numneigh[i];
rand_type rand_gen = rand_pool.get_state();
for (jj = 0; jj < jnum; jj++) {
j = d_neighbors(i,jj);
factor_dpd = special_lj[sbmask(j)];
j &= NEIGHMASK;
delx = xtmp - x(j,0);
dely = ytmp - x(j,1);
delz = ztmp - x(j,2);
rsq = delx*delx + dely*dely + delz*delz;
jtype = type(j);
if (rsq < d_cutsq(itype,jtype)) {
r = sqrt(rsq);
if (r < EPSILON) continue; // r can be 0.0 in DPD systems
rinv = 1.0/r;
delvx = vxtmp - v(j,0);
delvy = vytmp - v(j,1);
delvz = vztmp - v(j,2);
dot = delx*delvx + dely*delvy + delz*delvz;
wd = 1.0 - r/params(itype,jtype).cut;
randnum = rand_gen.normal();
// conservative force
fpair = params(itype,jtype).a0*wd;
// drag force - parallel
fpair -= params(itype,jtype).gamma*wd*wd*dot*rinv;
// random force - parallel
fpair += params(itype,jtype).sigma*wd*randnum*dtinvsqrt;
fpair *= factor_dpd*rinv;
fx += fpair*delx;
fy += fpair*dely;
fz += fpair*delz;
if ((neighflag==HALF || neighflag==HALFTHREAD) && (NEWTON_PAIR || j < nlocal) ) {
a_f(j,0) -= fpair*delx;
a_f(j,1) -= fpair*dely;
a_f(j,2) -= fpair*delz;
}
if (eflag) {
// unshifted eng of conservative term:
// evdwl = -a0[itype][jtype]*r * (1.0-0.5*r/cut[itype][jtype]);
// eng shifted to 0.0 at cutoff
evdwl = 0.5*params(itype,jtype).a0*params(itype,jtype).cut* wd*wd;
evdwl *= factor_dpd;
if (EVFLAG)
ev.evdwl += (((NEIGHFLAG==HALF || NEIGHFLAG==HALFTHREAD) && (NEWTON_PAIR||(j<nlocal)))?1.0:0.5)*evdwl;
}
if (EVFLAG)
this->template ev_tally<NEIGHFLAG,NEWTON_PAIR>(ev,i,j,evdwl,fpair,delx,dely,delz);
}
}
a_f(i,0) += fx;
a_f(i,1) += fy;
a_f(i,2) += fz;
rand_pool.free_state(rand_gen);
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
template<int NEIGHFLAG, int NEWTON_PAIR>
KOKKOS_INLINE_FUNCTION
void PairDPDKokkos<DeviceType>::ev_tally(EV_FLOAT &ev, const int &i, const int &j,
const F_FLOAT &epair, const F_FLOAT &fpair, const F_FLOAT &delx,
const F_FLOAT &dely, const F_FLOAT &delz) const
{
const int EFLAG = eflag;
const int VFLAG = vflag_either;
// The eatom and vatom arrays are atomic for Half/Thread neighbor style
Kokkos::View<E_FLOAT*, typename DAT::t_efloat_1d::array_layout,typename KKDevice<DeviceType>::value,Kokkos::MemoryTraits<AtomicF<NEIGHFLAG>::value> > a_eatom = k_eatom.view<DeviceType>();
Kokkos::View<F_FLOAT*[6], typename DAT::t_virial_array::array_layout,typename KKDevice<DeviceType>::value,Kokkos::MemoryTraits<AtomicF<NEIGHFLAG>::value> > a_vatom = k_vatom.view<DeviceType>();
if (EFLAG) {
if (eflag_atom) {
const E_FLOAT epairhalf = 0.5 * epair;
if (NEIGHFLAG!=FULL) {
if (NEWTON_PAIR || i < nlocal) a_eatom[i] += epairhalf;
if (NEWTON_PAIR || j < nlocal) a_eatom[j] += epairhalf;
} else {
a_eatom[i] += epairhalf;
}
}
}
if (VFLAG) {
const E_FLOAT v0 = delx*delx*fpair;
const E_FLOAT v1 = dely*dely*fpair;
const E_FLOAT v2 = delz*delz*fpair;
const E_FLOAT v3 = delx*dely*fpair;
const E_FLOAT v4 = delx*delz*fpair;
const E_FLOAT v5 = dely*delz*fpair;
if (vflag_global) {
if (NEIGHFLAG!=FULL) {
if (NEWTON_PAIR || i < nlocal) {
ev.v[0] += 0.5*v0;
ev.v[1] += 0.5*v1;
ev.v[2] += 0.5*v2;
ev.v[3] += 0.5*v3;
ev.v[4] += 0.5*v4;
ev.v[5] += 0.5*v5;
}
if (NEWTON_PAIR || j < nlocal) {
ev.v[0] += 0.5*v0;
ev.v[1] += 0.5*v1;
ev.v[2] += 0.5*v2;
ev.v[3] += 0.5*v3;
ev.v[4] += 0.5*v4;
ev.v[5] += 0.5*v5;
}
} else {
ev.v[0] += 0.5*v0;
ev.v[1] += 0.5*v1;
ev.v[2] += 0.5*v2;
ev.v[3] += 0.5*v3;
ev.v[4] += 0.5*v4;
ev.v[5] += 0.5*v5;
}
}
if (vflag_atom) {
if (NEIGHFLAG!=FULL) {
if (NEWTON_PAIR || i < nlocal) {
a_vatom(i,0) += 0.5*v0;
a_vatom(i,1) += 0.5*v1;
a_vatom(i,2) += 0.5*v2;
a_vatom(i,3) += 0.5*v3;
a_vatom(i,4) += 0.5*v4;
a_vatom(i,5) += 0.5*v5;
}
if (NEWTON_PAIR || j < nlocal) {
a_vatom(j,0) += 0.5*v0;
a_vatom(j,1) += 0.5*v1;
a_vatom(j,2) += 0.5*v2;
a_vatom(j,3) += 0.5*v3;
a_vatom(j,4) += 0.5*v4;
a_vatom(j,5) += 0.5*v5;
}
} else {
a_vatom(i,0) += 0.5*v0;
a_vatom(i,1) += 0.5*v1;
a_vatom(i,2) += 0.5*v2;
a_vatom(i,3) += 0.5*v3;
a_vatom(i,4) += 0.5*v4;
a_vatom(i,5) += 0.5*v5;
}
}
}
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
void PairDPDKokkos<DeviceType>::allocate()
{
PairDPD::allocate();
int n = atom->ntypes;
memory->destroy(cutsq);
memoryKK->create_kokkos(k_cutsq,cutsq,n+1,n+1,"pair:cutsq");
d_cutsq = k_cutsq.template view<DeviceType>();
k_params = Kokkos::DualView<params_dpd**,Kokkos::LayoutRight,DeviceType>("PairDPD::params",n+1,n+1);
params = k_params.template view<DeviceType>();
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
KOKKOS_INLINE_FUNCTION
int PairDPDKokkos<DeviceType>::sbmask(const int& j) const {
return j >> SBBITS & 3;
}
/* ----------------------------------------------------------------------
init for one type pair i,j and corresponding j,i
------------------------------------------------------------------------- */
template<class DeviceType>
double PairDPDKokkos<DeviceType>::init_one(int i, int j)
{
double cutone = PairDPD::init_one(i,j);
k_params.h_view(i,j).cut = cut[i][j];
k_params.h_view(i,j).a0 = a0[i][j];
k_params.h_view(i,j).gamma = gamma[i][j];
k_params.h_view(i,j).sigma = sigma[i][j];
k_params.h_view(j,i) = k_params.h_view(i,j);
k_params.template modify<LMPHostType>();
k_cutsq.h_view(i,j) = cutone*cutone;
k_cutsq.h_view(j,i) = k_cutsq.h_view(i,j);
k_cutsq.template modify<LMPHostType>();
return cutone;
}
namespace LAMMPS_NS {
template class PairDPDKokkos<LMPDeviceType>;
#ifdef LMP_KOKKOS_GPU
template class PairDPDKokkos<LMPHostType>;
#endif
}

View File

@ -0,0 +1,126 @@
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
https://www.lammps.org/, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
#ifdef PAIR_CLASS
// clang-format off
PairStyle(dpd/kk,PairDPDKokkos<LMPDeviceType>);
PairStyle(dpd/kk/device,PairDPDKokkos<LMPDeviceType>);
PairStyle(dpd/kk/host,PairDPDKokkos<LMPHostType>);
// clang-format on
#else
#ifndef LMP_PAIR_DPD_KOKKOS_H
#define LMP_PAIR_DPD_KOKKOS_H
#include "pair_dpd.h"
#include "pair_kokkos.h"
#include "kokkos_type.h"
#if !defined(DPD_USE_RAN_MARS) && !defined(DPD_USE_Random_XorShift64) && !defined(Random_XorShift1024)
#define DPD_USE_Random_XorShift64
#endif
#ifdef DPD_USE_RAN_MARS
#include "rand_pool_wrap_kokkos.h"
#else
#include "Kokkos_Random.hpp"
#endif
namespace LAMMPS_NS {
template<class DeviceType>
class PairDPDKokkos : public PairDPD {
public:
typedef DeviceType device_type;
typedef ArrayTypes<DeviceType> AT;
typedef EV_FLOAT value_type;
PairDPDKokkos(class LAMMPS*);
~PairDPDKokkos() override;
void allocate() override;
void init_style() override;
double init_one(int i, int j) override;
void compute(int, int) override;
struct params_dpd {
KOKKOS_INLINE_FUNCTION
params_dpd() {cut=a0=gamma=sigma=0;}
KOKKOS_INLINE_FUNCTION
params_dpd(int /*i*/) {cut=a0=gamma=sigma=0;}
F_FLOAT cutsq,cut,a0,gamma,sigma;
};
template<int NEIGHFLAG, int NEWTON_PAIR, int EVFLAG>
struct TagDPDKokkos{};
template<int NEIGHFLAG, int NEWTON_PAIR, int EVFLAG>
KOKKOS_INLINE_FUNCTION
void operator () (TagDPDKokkos<NEIGHFLAG,NEWTON_PAIR,EVFLAG>, const int &i) const;
template<int NEIGHFLAG, int NEWTON_PAIR, int EVFLAG>
KOKKOS_INLINE_FUNCTION
void operator () (TagDPDKokkos<NEIGHFLAG,NEWTON_PAIR,EVFLAG>, const int &i, EV_FLOAT&) const;
template<int NEIGHFLAG, int NEWTON_PAIR>
KOKKOS_INLINE_FUNCTION
void ev_tally(EV_FLOAT &ev, const int &i, const int &j,
const F_FLOAT &epair, const F_FLOAT &fpair, const F_FLOAT &delx,
const F_FLOAT &dely, const F_FLOAT &delz) const;
private:
double special_lj[4];
int eflag,vflag;
int neighflag, nlocal,newton_pair;
double dtinvsqrt;
#ifdef DPD_USE_RAN_MARS
RandPoolWrap rand_pool;
typedef RandWrap rand_type;
#elif defined(DPD_USE_Random_XorShift64)
Kokkos::Random_XorShift64_Pool<DeviceType> rand_pool;
typedef typename Kokkos::Random_XorShift64_Pool<DeviceType>::generator_type rand_type;
#elif defined(DPD_USE_Random_XorShift1024)
Kokkos::Random_XorShift1024_Pool<DeviceType> rand_pool;
typedef typename Kokkos::Random_XorShift1024_Pool<DeviceType>::generator_type rand_type;
#endif
typename AT::t_x_array_randomread x;
typename AT::t_x_array_randomread v;
typename AT::t_f_array f;
typename AT::t_int_1d_randomread type;
typename AT::t_neighbors_2d d_neighbors;
typename AT::t_int_1d_randomread d_ilist;
typename AT::t_int_1d_randomread d_numneigh;
typename AT::tdual_ffloat_2d k_cutsq;
typename AT::t_ffloat_2d d_cutsq;
Kokkos::DualView<params_dpd**,Kokkos::LayoutRight,DeviceType> k_params;
typename Kokkos::DualView<params_dpd**,
Kokkos::LayoutRight,DeviceType>::t_dev_const_um params;
DAT::tdual_efloat_1d k_eatom;
DAT::tdual_virial_array k_vatom;
typename AT::t_efloat_1d d_eatom;
typename AT::t_virial_array d_vatom;
KOKKOS_INLINE_FUNCTION
int sbmask(const int& j) const;
friend void pair_virial_fdotr_compute<PairDPDKokkos>(PairDPDKokkos*);
};
}
#endif
#endif