19 #ifndef _SF_ABSTRACT_VECTOR_H
20 #define _SF_ABSTRACT_VECTOR_H
53 template<
class T,
class S>
77 ltype inp_layout) = 0;
94 virtual void init(
int igsize,
int ilsize,
int idpn = 1,
ltype ilayout =
unset) = 0;
122 MPI_Allreduce(MPI_IN_PLACE, &N, 1, MPI_INT, MPI_SUM,
mesh->
comm);
142 return std::tuple<int, int>(N, n);
167 virtual void set(
const S val) = 0;
177 virtual void set(
const T idx,
const S val) = 0;
195 virtual S
get(
const T idx) = 0;
357 virtual S
mag()
const = 0;
364 virtual S
sum()
const = 0;
371 virtual S
min()
const = 0;
445 MPI_Comm comm = this->mesh != NULL ? this->mesh->
comm : PETSC_COMM_WORLD;
446 MPI_Comm_rank(comm, &rank);
447 MPI_Comm_size(comm, &size);
449 int glb_size = this->
gsize();
450 int loc_size = this->
lsize();
457 fd = fopen(file,
"w");
461 MPI_Allreduce(MPI_IN_PLACE, &err, 1, MPI_INT, MPI_MAX, comm);
471 fprintf(fd,
"%d\n", glb_size);
473 for(
int i=0; i<loc_size/this->
dpn; i++) {
474 for(
int j=0; j<this->
dpn; j++)
475 nwr += fprintf(fd,
"%f ", p[i*this->dpn+j]);
476 nwr += fprintf(fd,
"\n");
480 for(
int pid=1; pid < size; pid++)
485 MPI_Recv(&rsize, 1, MPI_INT, pid,
SF_MPITAG, comm, &stat);
487 MPI_Recv(wbuff.
data(), rsize*
sizeof(S), MPI_BYTE, pid,
SF_MPITAG, comm, &stat);
489 for(
int i=0; i<rsize/this->
dpn; i++) {
490 for(
int j=0; j<this->
dpn; j++)
491 nwr += fprintf(fd,
"%f ", wbuff[i*this->dpn+j]);
492 nwr += fprintf(fd,
"\n");
498 MPI_Send(&loc_size, 1, MPI_INT, 0,
SF_MPITAG, comm);
499 MPI_Send(p, loc_size*
sizeof(S), MPI_BYTE, 0,
SF_MPITAG, comm);
503 MPI_Bcast(&nwr, 1, MPI_LONG, 0, comm);
519 MPI_Comm comm = this->mesh != NULL ? this->mesh->
comm : PETSC_COMM_WORLD;
520 MPI_Comm_rank(comm, &rank); MPI_Comm_size(comm, &size);
522 long int loc_size = this->
lsize();
526 for(
long int i=0; i<loc_size; i++) buff[i] = p[i];
539 MPI_Comm comm = this->mesh != NULL ? this->mesh->
comm : PETSC_COMM_WORLD;
540 int rank; MPI_Comm_rank(comm, &rank);
546 fd = fopen(file.c_str(),
"w");
551 MPI_Allreduce(MPI_IN_PLACE, &error, 1, MPI_INT, MPI_SUM, comm);
554 nwr = this->write_binary<V>(fd);
567 MPI_Comm comm = this->mesh != NULL ? this->mesh->
comm : PETSC_COMM_WORLD;
569 size_t loc_size = this->
lsize();
575 for(
size_t i=0; i<loc_size; i++)
586 MPI_Comm comm =
mesh != NULL ?
mesh->
comm : PETSC_COMM_WORLD;
589 int rank; MPI_Comm_rank(comm, &rank);
595 fd = fopen(file.c_str(),
"r");
600 MPI_Allreduce(MPI_IN_PLACE, &error, 1, MPI_INT, MPI_SUM, comm);
603 nrd = read_binary<V>(fd);
615 MPI_Comm comm = this->mesh != NULL ? this->mesh->
comm : PETSC_COMM_WORLD;
617 size_t loc_size = this->
lsize();
626 MPI_Comm comm = this->mesh != NULL ? this->mesh->
comm : PETSC_COMM_WORLD;
627 int rank; MPI_Comm_rank(comm, &rank);
635 fd = fopen(file.c_str(),
"r");
640 MPI_Allreduce(MPI_IN_PLACE, &error, 1, MPI_INT, MPI_SUM, comm);
657 template<
class T,
class S>
699 template<
class T,
class S>
711 template<
class T,
class S>
723 template<
class T,
class S>
762 const size_t gsize_a,
763 const size_t gsize_b,
766 assert(_registry.
count(spec) == 0);
769 _registry[spec] = sc;
778 for(
size_t i=0; i<nbr_a.
size(); i++)
779 for(
short j=0; j<dpn; j++) idx_a[i*dpn+j] = nbr_a[i]*dpn+j;
781 for(
size_t i=0; i<nbr_b.
size(); i++)
782 for(
short j=0; j<dpn; j++) idx_b[i*dpn+j] = nbr_b[i]*dpn+j;
784 err += ISCreateGeneral(PETSC_COMM_WORLD, idx_a.
size(), idx_a.
data(), PETSC_COPY_VALUES, &is_a);
785 err += ISCreateGeneral(PETSC_COMM_WORLD, idx_b.
size(), idx_b.
data(), PETSC_COPY_VALUES, &is_b);
787 err += ISSetPermutation(is_b);
790 err += VecCreateMPI(PETSC_COMM_WORLD, idx_a.
size(), gsize_a*dpn, &a);
791 err += VecCreateMPI(PETSC_COMM_WORLD, idx_b.
size(), gsize_b*dpn, &sc->
b_buff);
793 err += VecScatterCreate(a, is_a, sc->
b_buff, is_b, &sc->
vec_sc);
795 err += VecDestroy(&a);
796 err += ISDestroy (&is_a);
797 err += ISDestroy (&is_b);
800 fprintf(stderr,
"%s : an error has occurred! Scattering spec: Mesh ID %d, permutation ID %d, numbering %d, dpn %d!\n",
801 __func__,
int(spec.
v1),
int(spec.
v2),
int(spec.
v3),
int(spec.
v4));
825 assert(_registry.
count(spec) == 0);
828 assert(idx_a.
size() == idx_b.
size());
831 _registry[spec] = sc;
839 T lsize_a = layout_a[rank+1] - layout_a[rank];
840 T lsize_b = layout_b[rank+1] - layout_b[rank];
841 T gsize_a = layout_a[layout_a.
size()-1];
842 T gsize_b = layout_b[layout_a.
size()-1];
846 for(
size_t i=0; i<idx_a.
size(); i++)
847 for(
short j=0; j<dpn; j++) sidx_a[i*dpn+j] = idx_a[i]*dpn+j;
849 for(
size_t i=0; i<idx_b.
size(); i++)
850 for(
short j=0; j<dpn; j++) sidx_b[i*dpn+j] = idx_b[i]*dpn+j;
852 err += ISCreateGeneral(PETSC_COMM_WORLD, sidx_a.
size(), sidx_a.
data(), PETSC_COPY_VALUES, &is_a);
853 err += ISCreateGeneral(PETSC_COMM_WORLD, sidx_b.
size(), sidx_b.
data(), PETSC_COPY_VALUES, &is_b);
855 err += VecCreateMPI(PETSC_COMM_WORLD, lsize_a*dpn, gsize_a*dpn, &a);
856 err += VecCreateMPI(PETSC_COMM_WORLD, lsize_b*dpn, gsize_b*dpn, &sc->
b_buff);
858 err += VecScatterCreate(a, is_a, sc->
b_buff, is_b, &sc->
vec_sc);
860 err += VecDestroy(&a);
861 err += ISDestroy (&is_a);
862 err += ISDestroy (&is_b);
865 fprintf(stderr,
"%s : an error has occurred! Scattering spec: To ID %d, From ID %d, numbering %d, dpn %d!\n",
866 __func__,
int(spec.
v1),
int(spec.
v2),
int(spec.
v3),
int(spec.
v4));
889 if(_registry.
count(spec))
890 ret = _registry[spec];
901 for(it = _registry.
begin(); it != _registry.
end(); ++it) {
#define SF_MPITAG
the MPI tag when communicating
Classes and algorithms related to the layout of distributed meshes.
Simple utility functions for parallel data.
size_t read_ascii(FILE *fd)
size_t read_binary(std::string file)
virtual void get(const vector< T > &idx, S *out)=0
virtual void apply_scattering(scattering &sc, bool fwd)=0
virtual void operator-=(const abstract_vector< T, S > &vec)=0
virtual void release_ptr(S *&p)=0
virtual std::string to_string() const =0
virtual void operator*=(const abstract_vector< T, S > &vec)=0
virtual void init(int igsize, int ilsize, int idpn=1, ltype ilayout=unset)=0
virtual T gsize() const =0
size_t write_ascii(const char *file, bool write_header)
std::tuple< int, int > init_common(const meshdata< mesh_int_t, mesh_real_t > &imesh, int idpn, ltype inp_layout)
size_t read_binary(FILE *fd)
virtual void operator/=(const S sca)=0
virtual void forward(abstract_vector< T, S > &out, scattering &sc, bool add=false)=0
virtual void set(const vector< T > &idx, const S val)=0
virtual bool is_init() const =0
ltype layout
used vector layout (nodal, algebraic, unset)
virtual void set(const vector< T > &idx, const vector< S > &vals, const bool additive=false)=0
virtual bool equals(const abstract_vector< T, S > &rhs) const =0
virtual void const_release_ptr(const S *&p) const =0
size_t read_ascii(std::string file)
virtual void operator=(const abstract_vector< T, S > &rhs)=0
virtual void deep_copy(const abstract_vector< T, S > &v)=0
size_t write_binary(std::string file)
write binary. Open file descriptor myself.
virtual void get_ownership_range(T &start, T &stop) const =0
virtual S get(const T idx)=0
virtual void shallow_copy(const abstract_vector< T, S > &v)=0
virtual void add_scaled(const abstract_vector< T, S > &vec, S k)=0
virtual void overshadow(const abstract_vector< T, S > &sub, bool member, int offset, int sz, bool share)=0
virtual void set(const S val)=0
virtual ~abstract_vector()=default
virtual void operator+=(const abstract_vector< T, S > &vec)=0
virtual void operator=(const vector< S > &rhs)=0
virtual void set(const T idx, const S val)=0
size_t write_binary(FILE *fd)
Write a vector to HD in binary. File descriptor is already set up.
virtual void backward(abstract_vector< T, S > &out, scattering &sc, bool add=false)=0
virtual T lsize() const =0
virtual void init(const abstract_vector< T, S > &vec)=0
int dpn
d.o.f. per mesh vertex.
virtual const S * const_ptr() const =0
virtual void operator*=(const S sca)=0
virtual void finish_assembly()=0
virtual void init(const meshdata< mesh_int_t, mesh_real_t > &imesh, int idpn, ltype inp_layout)=0
const meshdata< mesh_int_t, mesh_real_t > * mesh
the connected mesh
virtual void operator+=(S c)=0
virtual S dot(const abstract_vector< T, S > &v) const =0
overlapping_layout< T > pl
nodal parallel layout
size_t l_numelem
local number of elements
size_t l_numpts
local number of points
size_t g_numelem
global number of elements
MPI_Comm comm
the parallel mesh is defined on a MPI world
The scatterer registry class.
scattering * register_scattering(const quadruple< int > spec, const vector< T > &layout_a, const vector< T > &layout_b, const vector< T > &idx_a, const vector< T > &idx_b, const int rank, const int dpn)
Register a scattering.
scattering * register_permutation(const quadruple< int > spec, const vector< T > &nbr_a, const vector< T > &nbr_b, const size_t gsize_a, const size_t gsize_b, const short dpn)
Register a permutation scattering.
scattering * get_scattering(const quadruple< int > spec)
Access an previously registered scattering.
void free_scatterings()
Free the registered scatterings.
Container for a PETSc VecScatter.
Vec b_buff
A buffer vector which also defines the parallel layout of the "b" side.
void operator()(abstract_vector< T, S > &v, bool fwd)
Apply the scattering on a data vector.
void forward(abstract_vector< T, S > &in, abstract_vector< T, S > &out, bool add=false)
Forward scattering.
VecScatter vec_sc
The scatterer.
void backward(abstract_vector< T, S > &in, abstract_vector< T, S > &out, bool add=false)
Backward scattering.
A vector storing arbitrary data.
size_t size() const
The current size of the vector.
void resize(size_t n)
Resize a vector.
T * data()
Pointer to the vector's start.
hm_int count(const K &key) const
Check if key exists.
Classes similar to unordered_set and unordered_map, but with better performance.
void treat_file_open_error(const char *file, const char *caller, const int errnum, const bool do_exit, int rank)
treat a file open error by displaying the errnum string interpretation and the caller
bool is_init(const abstract_vector< T, S > *v)
size_t root_read(FILE *fd, vector< V > &vec, MPI_Comm comm)
Read binary data into a vector.
size_t root_read_ascii(FILE *fd, vector< V > &vec, MPI_Comm comm, bool int_data)
Read binary data into a vector.
size_t root_write(FILE *fd, const vector< V > &vec, MPI_Comm comm)
Write vector data binary to disk.