void cow_dfield_sampleexecute(cow_dfield *f) // ----------------------------------------------------------------------------- // Samples `N` points on the global domain. `N` may vary between processes, or // be equal to zero. Points are output permuted relative to the input // coordinates, so take care to use `xout` as the coordinates associated with // samples in `P`. // // f: IN data field instance to sample // x: IN list of input coordinates at which to sample f's data (N x 3) // N: IN number of points to sample // xout: OUT locations of returned samples, permutation of x (N x 3) // P: OUT list of filled samples (N x Q) where Q = f->n_members // ----------------------------------------------------------------------------- { double *xout = (double*) malloc(f->samplecoordslen * 3 * sizeof(double)); double *xin = f->samplecoords; int N = f->samplecoordslen; double *P = f->sampleresult; int mode = f->samplemode; if (cow_mpirunning()) { _rem(f, xin, N, xout, P, mode); } else { _loc(f, xin, N, xout, P, mode); } memcpy(f->samplecoords, xout, f->samplecoordslen * 3 * sizeof(double)); free(xout); }
void _io_domain_commit(cow_domain *d) { #if (COW_HDF5) for (int n=0; n<d->n_dims; ++n) { d->L_nint_h5[n] = d->L_nint[n]; // Selection size, target and destination d->L_ntot_h5[n] = d->L_ntot[n]; // Memory space total size d->L_strt_h5[n] = d->L_strt[n]; // Memory space selection start d->G_ntot_h5[n] = d->G_ntot[n]; // Global space total size d->G_strt_h5[n] = d->G_strt[n]; // Global space selection start } // Here we create the following property lists: // // file access property list ........ for the call to H5Fopen // dset creation property list ........ for the call to H5Dcreate // dset transfer property list ........ for the call to H5Dwrite // --------------------------------------------------------------------------- d->fapl = H5Pcreate(H5P_FILE_ACCESS); d->dcpl = H5Pcreate(H5P_DATASET_CREATE); d->dxpl = H5Pcreate(H5P_DATASET_XFER); #if (COW_HDF5_MPI && COW_MPI) if (cow_mpirunning()) { H5Pset_fapl_mpio(d->fapl, d->mpi_cart, MPI_INFO_NULL); } #endif // COW_HDF5_MPI && COW_MPI #endif // COW_HDF5 }
void cow_histogram_seal(cow_histogram *h) { if (!h->committed || h->sealed) return; #if (COW_MPI) if (cow_mpirunning()) { int nbins = h->nbinsx * h->nbinsy; MPI_Comm c = h->comm; MPI_Allreduce(MPI_IN_PLACE, h->weight, nbins, MPI_DOUBLE, MPI_SUM, c); MPI_Allreduce(MPI_IN_PLACE, h->counts, nbins, MPI_LONG, MPI_SUM, c); MPI_Allreduce(MPI_IN_PLACE, &h->totcounts, 1, MPI_LONG, MPI_SUM, c); } #endif h->sealed = 1; _filloutput(h); }
void cow_domain_setcollective(cow_domain *d, int mode) { #if (COW_HDF5 && COW_HDF5_MPI) if (mode && !cow_mpirunning()) { printf("[%s] requested collective without MPI running: " "revert to independent\n", MODULE); mode = 0; } if (mode) { printf("[%s] setting HDF5 io mode to collective\n", MODULE); H5Pset_dxpl_mpio(d->dxpl, H5FD_MPIO_COLLECTIVE); } else { printf("[%s] setting HDF5 io mode to independent\n", MODULE); H5Pset_dxpl_mpio(d->dxpl, H5FD_MPIO_INDEPENDENT); } #endif }
void cow_histogram_dumpascii(cow_histogram *h, char *fn) // ----------------------------------------------------------------------------- // Dumps the histogram as ascii to the file named `fn`. Synchronizes it across // processes before doing so. The function uses rank 0 to do the write. // ----------------------------------------------------------------------------- { if (!(h->committed && h->sealed)) { return; } #if (COW_MPI) if (cow_mpirunning()) { int rank; MPI_Comm_rank(h->comm, &rank); if (rank != 0) { return; } } #endif FILE *file = fopen(fn, "w"); if (file == NULL) { printf("[%s] could not open file %s\n", __FILE__, fn); return; } else { printf("[%s] writing histogram as ASCII table to %s\n", MODULE, fn); } if (h->n_dims == 1) { for (int n=0; n<h->nbinsx; ++n) { fprintf(file, "%f %f\n", h->binlocx[n], h->binvalv[n]); } } else if (h->n_dims == 2) { for (int nx=0; nx<h->nbinsx; ++nx) { for (int ny=0; ny<h->nbinsy; ++ny) { fprintf(file, "%f %f %f\n", h->binlocx[nx], h->binlocy[ny], h->binvalv[nx * h->nbinsy + ny]); } } } fclose(file); }
FFT_DATA *_fwd(cow_dfield *f, double *fx, int start, int stride) { FFT_DATA *Fk = NULL; FFT_DATA *Fx = NULL; if (cow_mpirunning()) { #if (COW_MPI) int nbuf; long long ntot = cow_domain_getnumglobalzones(f->domain, COW_ALL_DIMS); struct fft_plan_3d *plan = call_fft_plan_3d(f->domain, &nbuf); Fx = (FFT_DATA*) malloc(nbuf * sizeof(FFT_DATA)); Fk = (FFT_DATA*) malloc(nbuf * sizeof(FFT_DATA)); for (int n=0; n<nbuf; ++n) { Fx[n][0] = fx[stride * n + start] / ntot; Fx[n][1] = 0.0; } fft_3d(Fx, Fk, FFT_FWD, plan); free(Fx); fft_3d_destroy_plan(plan); #endif // COW_MPI } else { int nbuf = cow_domain_getnumlocalzonesinterior(f->domain, COW_ALL_DIMS); long long ntot = cow_domain_getnumglobalzones(f->domain, COW_ALL_DIMS); Fx = (FFT_DATA*) malloc(nbuf * sizeof(FFT_DATA)); Fk = (FFT_DATA*) malloc(nbuf * sizeof(FFT_DATA)); for (int n=0; n<nbuf; ++n) { Fx[n][0] = fx[stride * n + start] / ntot; Fx[n][1] = 0.0; } int *N = f->domain->L_nint; fftw_plan plan = fftw_plan_many_dft(3, N, 1, Fx, NULL, 1, 0, Fk, NULL, 1, 0, FFTW_FORWARD, FFTW_ESTIMATE); fftw_execute(plan); fftw_destroy_plan(plan); free(Fx); } return Fk; }
double *_rev(cow_domain *d, FFT_DATA *Fk) { FFT_DATA *Fx = NULL; double *fx = NULL; if (cow_mpirunning()) { #if (COW_MPI) int nbuf; long long ntot = cow_domain_getnumglobalzones(d, COW_ALL_DIMS); struct fft_plan_3d *plan = call_fft_plan_3d(d, &nbuf); fx = (double*) malloc(nbuf * sizeof(double)); Fx = (FFT_DATA*) malloc(nbuf * sizeof(FFT_DATA)); fft_3d(Fk, Fx, FFT_REV, plan); for (int n=0; n<nbuf; ++n) { fx[n] = Fx[n][0] / ntot; } free(Fx); fft_3d_destroy_plan(plan); #endif // COW_MPI } else { int nbuf = cow_domain_getnumlocalzonesinterior(d, COW_ALL_DIMS); long long ntot = cow_domain_getnumglobalzones(d, COW_ALL_DIMS); fx = (double*) malloc(nbuf * sizeof(double)); Fx = (FFT_DATA*) malloc(nbuf * sizeof(FFT_DATA)); int *N = d->L_nint; fftw_plan plan = fftw_plan_many_dft(3, N, 1, Fk, NULL, 1, 0, Fx, NULL, 1, 0, FFTW_BACKWARD, FFTW_ESTIMATE); fftw_execute(plan); for (int n=0; n<nbuf; ++n) { fx[n] = Fx[n][0] / ntot; } free(Fx); fftw_destroy_plan(plan); } return fx; }
void _io_read(cow_dfield *f, char *fname) { #if (COW_HDF5) cow_domain *d = f->domain; char **pnames = f->members; void *data = f->data; char *gname = f->name; int n_memb = f->n_members; int n_dims = d->n_dims; hsize_t *L_nint = d->L_nint_h5; hsize_t *G_strt = d->G_strt_h5; hsize_t *G_ntot = d->G_ntot_h5; hsize_t ndp1 = n_dims + 1; hsize_t l_nint[4]; hsize_t l_ntot[4]; hsize_t l_strt[4]; hsize_t stride[4]; for (int i=0; i<n_dims; ++i) { l_nint[i] = d->L_nint[i]; // Selection size, target and destination l_ntot[i] = d->L_ntot[i]; // Memory space total size l_strt[i] = d->L_strt[i]; // Memory space selection start stride[i] = 1; } l_nint[ndp1 - 1] = 1; l_ntot[ndp1 - 1] = n_memb; stride[ndp1 - 1] = n_memb; // The loop over processors is needed if COW_MPI support is enabled and // COW_HDF5_MPI is not. If either COW_MPI is disabled, or COW_HDF5_MPI is // enabled, then the write calls occur without the loop. // --------------------------------------------------------------------------- #if (!COW_HDF5_MPI && COW_MPI) for (int rank=0; rank<d->cart_size; ++rank) { if (rank == d->cart_rank) { #endif hid_t file = H5Fopen(fname, H5F_ACC_RDONLY, d->fapl); hid_t memb = H5Gopen(file, gname, H5P_DEFAULT); hid_t mspc = H5Screate_simple(ndp1, l_ntot, NULL); hid_t fspc = H5Screate_simple(n_dims, G_ntot, NULL); for (int n=0; n<n_memb; ++n) { hid_t dset = H5Dopen(memb, pnames[n], H5P_DEFAULT); l_strt[ndp1 - 1] = n; H5Sselect_hyperslab(mspc, H5S_SELECT_SET, l_strt, stride, l_nint, NULL); H5Sselect_hyperslab(fspc, H5S_SELECT_SET, G_strt, NULL, L_nint, NULL); H5Dread(dset, H5T_NATIVE_DOUBLE, mspc, fspc, d->dxpl, data); H5Dclose(dset); } H5Sclose(fspc); H5Sclose(mspc); H5Gclose(memb); H5Fclose(file); #if (!COW_HDF5_MPI && COW_MPI) } if (cow_mpirunning()) { MPI_Barrier(d->mpi_cart); } } #endif // !COW_HDF5_MPI && COW_MPI #endif }
void _io_write(cow_dfield *f, char *fname) // ----------------------------------------------------------------------------- // This function uses a collective MPI-IO procedure to write the contents of // 'data' to the HDF5 file named 'fname', which is assumed to have been created // already. The dataset with name 'dname', which is being written to, must not // exist already. Chunking is enabled as per the ChunkSize variable, and is // disabled by default. Recommended chunk size is local subdomain size. This // will result in optimized read/write on the same decomposition layout, but // poor performance for different access patterns, for example the slabs used by // cluster-FFT functions. // // WARNING! // // All processors must define the same chunk size, the behavior of this function // is not defined otherwise. This implies that chunking should be disabled when // running on a strange number of cores, and subdomain sizes are non-uniform. // ----------------------------------------------------------------------------- { #if (COW_HDF5) cow_domain *d = f->domain; char **pnames = f->members; void *data = f->data; char *gname = f->name; int n_memb = f->n_members; int n_dims = d->n_dims; hsize_t *L_nint = d->L_nint_h5; hsize_t *G_strt = d->G_strt_h5; hsize_t *G_ntot = d->G_ntot_h5; hsize_t ndp1 = n_dims + 1; hsize_t l_nint[4]; hsize_t l_ntot[4]; hsize_t l_strt[4]; hsize_t stride[4]; for (int i=0; i<n_dims; ++i) { l_nint[i] = d->L_nint[i]; // Selection size, target and destination l_ntot[i] = d->L_ntot[i]; // Memory space total size l_strt[i] = d->L_strt[i]; // Memory space selection start stride[i] = 1; } l_nint[ndp1 - 1] = 1; l_ntot[ndp1 - 1] = n_memb; stride[ndp1 - 1] = n_memb; // The loop over processors is needed if COW_MPI support is enabled and // COW_HDF5_MPI is not. If either COW_MPI is disabled, or COW_HDF5_MPI is // enabled, then the write calls occur without the loop. // --------------------------------------------------------------------------- #if (!COW_HDF5_MPI && COW_MPI) for (int rank=0; rank<d->cart_size; ++rank) { if (rank == d->cart_rank) { #endif hid_t file = H5Fopen(fname, H5F_ACC_RDWR, d->fapl); hid_t memb = H5Gopen(file, gname, H5P_DEFAULT); hid_t mspc = H5Screate_simple(ndp1, l_ntot, NULL); hid_t fspc = H5Screate_simple(n_dims, G_ntot, NULL); for (int n=0; n<n_memb; ++n) { int dset = H5Dcreate(memb, pnames[n], H5T_NATIVE_DOUBLE, fspc, H5P_DEFAULT, d->dcpl, H5P_DEFAULT); l_strt[ndp1 - 1] = n; H5Sselect_hyperslab(mspc, H5S_SELECT_SET, l_strt, stride, l_nint, NULL); H5Sselect_hyperslab(fspc, H5S_SELECT_SET, G_strt, NULL, L_nint, NULL); H5Dwrite(dset, H5T_NATIVE_DOUBLE, mspc, fspc, d->dxpl, data); H5Dclose(dset); } H5Sclose(fspc); H5Sclose(mspc); H5Gclose(memb); H5Fclose(file); #if (!COW_HDF5_MPI && COW_MPI) } if (cow_mpirunning()) { MPI_Barrier(d->mpi_cart); } } #endif // !COW_HDF5_MPI && COW_MPI #endif }
void cow_histogram_dumphdf5(cow_histogram *h, char *fn, char *gn) // ----------------------------------------------------------------------------- // Dumps the histogram to the HDF5 file named `fn`, under the group // `gn`/h->nickname, gn may be NULL. The function uses rank 0 to do the write. // ----------------------------------------------------------------------------- { #if (COW_HDF5) if (!(h->committed && h->sealed)) { return; } char gname[1024]; int rank = 0; if (gn) { snprintf(gname, 1024, "%s/%s", gn, h->nickname); } else { snprintf(gname, 1024, "%s", h->nickname); } #if (COW_MPI) if (cow_mpirunning()) { MPI_Comm_rank(h->comm, &rank); } #endif if (rank == 0) { // ------------------------------------------------------------------------- // The write functions assume the file is already created. Have master // create the file if it's not there already. // ------------------------------------------------------------------------- FILE *testf = fopen(fn, "r"); hid_t fid; if (testf == NULL) { fid = H5Fcreate(fn, H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT); } else { fclose(testf); fid = H5Fopen(fn, H5F_ACC_RDWR, H5P_DEFAULT); } if (H5Lexists_safe(fid, gname)) { printf("[%s] writing histogram as HDF5 to %s/%s (clobber existing)\n", MODULE, fn, gname); H5Gunlink(fid, gname); } else { printf("[%s] writing histogram as HDF5 to %s/%s\n", MODULE, fn, gname); } hid_t gcpl = H5Pcreate(H5P_LINK_CREATE); H5Pset_create_intermediate_group(gcpl, 1); hid_t memb = H5Gcreate(fid, gname, gcpl, H5P_DEFAULT, H5P_DEFAULT); H5Pclose(gcpl); H5Gclose(memb); H5Fclose(fid); } else { return; } // Create a group to represent this histogram, and an attribute to name it // --------------------------------------------------------------------------- hid_t fid = H5Fopen(fn, H5F_ACC_RDWR, H5P_DEFAULT); hid_t grp = H5Gopen(fid, gname, H5P_DEFAULT); if (h->fullname != NULL) { hid_t aspc = H5Screate(H5S_SCALAR); hid_t strn = H5Tcopy(H5T_C_S1); H5Tset_size(strn, strlen(h->fullname)); hid_t attr = H5Acreate(grp, "fullname", strn, aspc, H5P_DEFAULT, H5P_DEFAULT); H5Awrite(attr, strn, h->fullname); // write the full name H5Aclose(attr); H5Tclose(strn); H5Sclose(aspc); } // Create the data sets in the group: binloc (bin centers) and binval (values) // --------------------------------------------------------------------------- double *binlocX = h->binlocx; double *binlocY = h->binlocy; double *binvalV = h->binvalv; hsize_t sizeX[1] = { h->nbinsx }; hsize_t sizeY[1] = { h->nbinsy }; hsize_t SizeX[1] = { h->nbinsx + 1 }; // to hold bin edges hsize_t SizeY[1] = { h->nbinsy + 1 }; // below, cap S/F refers to bin edges hsize_t sizeZ[2] = { h->nbinsx, h->nbinsy }; hid_t fspcZ = H5Screate_simple(h->n_dims, sizeZ, NULL); if (h->n_dims >= 1) { hid_t fspcX = H5Screate_simple(1, sizeX, NULL); hid_t FspcX = H5Screate_simple(1, SizeX, NULL); hid_t dsetbinX = H5Dcreate(grp, "binlocX", H5T_NATIVE_DOUBLE, fspcX, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); hid_t dsetedgX = H5Dcreate(grp, "binedgX", H5T_NATIVE_DOUBLE, FspcX, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); H5Dwrite(dsetbinX, H5T_NATIVE_DOUBLE, fspcX, fspcX, H5P_DEFAULT, binlocX); H5Dwrite(dsetedgX, H5T_NATIVE_DOUBLE, FspcX, FspcX, H5P_DEFAULT, h->bedgesx); H5Dclose(dsetbinX); H5Sclose(FspcX); H5Sclose(fspcX); } if (h->n_dims >= 2) { hid_t fspcY = H5Screate_simple(1, sizeY, NULL); hid_t FspcY = H5Screate_simple(1, SizeY, NULL); hid_t dsetbinY = H5Dcreate(grp, "binlocY", H5T_NATIVE_DOUBLE, fspcY, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); hid_t dsetedgY = H5Dcreate(grp, "binedgY", H5T_NATIVE_DOUBLE, FspcY, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); H5Dwrite(dsetbinY, H5T_NATIVE_DOUBLE, fspcY, fspcY, H5P_DEFAULT, binlocY); H5Dwrite(dsetedgY, H5T_NATIVE_DOUBLE, FspcY, FspcY, H5P_DEFAULT, h->bedgesy); H5Dclose(dsetbinY); H5Sclose(FspcY); H5Sclose(fspcY); } hid_t dsetvalV = H5Dcreate(grp, "binval", H5T_NATIVE_DOUBLE, fspcZ, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); H5Dwrite(dsetvalV, H5T_NATIVE_DOUBLE, fspcZ, fspcZ, H5P_DEFAULT, binvalV); H5Dclose(dsetvalV); H5Sclose(fspcZ); H5Gclose(grp); H5Fclose(fid); #endif }
cow_histogram *cow_histogram_new() { cow_histogram *h = (cow_histogram*) malloc(sizeof(cow_histogram)); cow_histogram hist = { .nbinsx = 1, .nbinsy = 1, .x0 = 0.0, .x1 = 1.0, .y0 = 0.0, .y1 = 1.0, .bedgesx = NULL, .bedgesy = NULL, .weight = NULL, .totcounts = 0, .counts = NULL, .nickname = NULL, .fullname = NULL, .binmode = COW_HIST_BINMODE_COUNTS, .spacing = COW_HIST_SPACING_LINEAR, .n_dims = 0, .committed = 0, .sealed = 0, .transform = NULL, .binlocx = NULL, .binlocy = NULL, .binvalv = NULL, #if (COW_MPI) .comm = MPI_COMM_WORLD, #endif } ; *h = hist; cow_histogram_setnickname(h, "histogram"); return h; } void cow_histogram_commit(cow_histogram *h) { if (h->committed) return; h->n_dims = h->nbinsy > 1 ? 2 : 1; if (h->n_dims == 1) { double dx = (h->x1 - h->x0) / h->nbinsx; h->bedgesx = (double*) malloc((h->nbinsx+1)*sizeof(double)); h->weight = (double*) malloc((h->nbinsx)*sizeof(double)); h->counts = (long*) malloc((h->nbinsx)*sizeof(long)); for (int n=0; n<h->nbinsx+1; ++n) { if (h->spacing == COW_HIST_SPACING_LOG) { h->bedgesx[n] = h->x0 * pow(h->x1 / h->x0, (double)n / h->nbinsx); } else if (h->spacing == COW_HIST_SPACING_LINEAR) { h->bedgesx[n] = h->x0 + n * dx; } } for (int n=0; n<h->nbinsx; ++n) { h->counts[n] = 0; h->weight[n] = 0.0; } } else if (h->n_dims == 2) { int nbins = h->nbinsx * h->nbinsy; double dx = (h->y1 - h->y0) / h->nbinsx; double dy = (h->y1 - h->y0) / h->nbinsy; h->bedgesx = (double*) malloc((h->nbinsx+1)*sizeof(double)); h->bedgesy = (double*) malloc((h->nbinsy+1)*sizeof(double)); h->weight = (double*) malloc(nbins*sizeof(double)); h->counts = (long*) malloc(nbins*sizeof(long)); for (int n=0; n<h->nbinsx+1; ++n) { if (h->spacing == COW_HIST_SPACING_LOG) { h->bedgesx[n] = h->x0 * pow(h->x1/h->x0, (double)n / h->nbinsx); } else if (h->spacing == COW_HIST_SPACING_LINEAR) { h->bedgesx[n] = h->x0 + n * dx; } } for (int n=0; n<h->nbinsy+1; ++n) { if (h->spacing == COW_HIST_SPACING_LOG) { h->bedgesy[n] = h->y0 * pow(h->y1/h->y0, (double)n / h->nbinsy); } else if (h->spacing == COW_HIST_SPACING_LINEAR) { h->bedgesy[n] = h->y0 + n * dy; } } for (int n=0; n<nbins; ++n) { h->counts[n] = 0; h->weight[n] = 0.0; } } #if (COW_MPI) if (cow_mpirunning()) { MPI_Comm_dup(h->comm, &h->comm); } #endif h->committed = 1; } void cow_histogram_del(cow_histogram *h) { #if (COW_MPI) if (h->committed && cow_mpirunning()) { MPI_Comm_free(&h->comm); } #endif free(h->bedgesx); free(h->bedgesy); free(h->weight); free(h->counts); free(h->nickname); free(h->fullname); free(h->binlocx); free(h->binlocy); free(h->binvalv); free(h); } void cow_histogram_setdomaincomm(cow_histogram *h, cow_domain *d) { #if (COW_MPI) if (h->committed) return; h->comm = d->mpi_cart; #endif } void cow_histogram_setbinmode(cow_histogram *h, int binmode) { if (h->committed || h->sealed) return; switch (binmode) { case COW_HIST_BINMODE_DENSITY: h->binmode = binmode; break; case COW_HIST_BINMODE_AVERAGE: h->binmode = binmode; break; case COW_HIST_BINMODE_COUNTS: h->binmode = binmode; break; default: printf("[%s] error: no such bin mode\n", MODULE); break; } } void cow_histogram_setspacing(cow_histogram *h, int spacing) { if (h->committed || h->sealed) return; switch (spacing) { case COW_HIST_SPACING_LINEAR: h->spacing = spacing; break; case COW_HIST_SPACING_LOG: h->spacing = spacing; break; default: printf("[%s] error: no such spacing\n", MODULE); break; } } void cow_histogram_setnbins(cow_histogram *h, int dim, int nbins) { if (h->committed || h->sealed) return; switch (dim) { case 0: h->nbinsx = nbins; break; case 1: h->nbinsy = nbins; break; case COW_ALL_DIMS: h->nbinsx = h->nbinsy = nbins; break; default: break; } } void cow_histogram_setlower(cow_histogram *h, int dim, double v0) { if (h->committed || h->sealed) return; switch (dim) { case 0: h->x0 = v0; break; case 1: h->y0 = v0; break; case COW_ALL_DIMS: h->x0 = h->y0 = v0; break; default: break; } } void cow_histogram_setupper(cow_histogram *h, int dim, double v1) { if (h->committed || h->sealed) return; switch (dim) { case 0: h->x1 = v1; break; case 1: h->y1 = v1; break; case COW_ALL_DIMS: h->x1 = h->y1 = v1; break; default: break; } } void cow_histogram_setfullname(cow_histogram *h, char *fullname) { h->fullname = (char*) realloc(h->fullname, strlen(fullname)+1); strcpy(h->fullname, fullname); }