#include #include #include #include #include #include #include #include #include #include #include #include "cctk.h" #include "cctk_Parameters.h" #include "AMRwriter.hh" #include "AmrGridReader.hh" #ifdef HDF4 # include "HDFIO.hh" #endif #ifdef HDF5 # include "H5IO.hh" #endif #include "IEEEIO.hh" #include "IO.hh" // Hack to stop FlexIO data type clash with LAM MPI #undef BYTE #undef CHAR #include "CactusBase/IOUtil/src/ioGH.h" #include "bbox.hh" #include "data.hh" #include "gdata.hh" #include "ggf.hh" #include "vect.hh" #include "carpet.hh" #include "ioflexio.hh" extern "C" { static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc,v 1.32 2003/09/23 12:32:12 cvs_anon Exp $"; CCTK_FILEVERSION(Carpet_CarpetIOFlexIO_ioflexio_cc); } namespace CarpetIOFlexIO { using namespace std; using namespace Carpet; // Variable definitions int GHExtension; int IOMethod; vector do_truncate; vector > last_output; static const char* GetStringParameter (const char* const parametername, const char* const fallback); static int GetIntParameter (const char* const parametername, int fallback); static bool CheckForVariable (const cGH* const cgh, const char* const varlist, const int vindex); static void SetFlag (int index, const char* optstring, void* arg); int CarpetIOFlexIOStartup () { CCTK_RegisterBanner ("AMR 3D FlexIO I/O provided by CarpetIOFlexIO"); GHExtension = CCTK_RegisterGHExtension("CarpetIOFlexIO"); CCTK_RegisterGHExtensionSetupGH (GHExtension, SetupGH); IOMethod = CCTK_RegisterIOMethod ("IOFlexIO"); CCTK_RegisterIOMethodOutputGH (IOMethod, OutputGH); CCTK_RegisterIOMethodOutputVarAs (IOMethod, OutputVarAs); CCTK_RegisterIOMethodTimeToOutput (IOMethod, TimeToOutput); CCTK_RegisterIOMethodTriggerOutput (IOMethod, TriggerOutput); return 0; } void* SetupGH (tFleshConfig* const fc, const int convLevel, cGH* const cgh) { DECLARE_CCTK_PARAMETERS; // Truncate all files if this is not a restart do_truncate.resize(CCTK_NumVars(), true); // No iterations have yet been output last_output.resize(maxreflevels); for (int rl=0; rl=0 && n=0 && group<(int)Carpet::arrdata.size()); const int n0 = CCTK_FirstVarIndexI(group); assert (n0>=0 && n0=0 && varrecovered || stat(filename, &fileinfo)!=0) { writer = 0; if (CCTK_Equals(out3D_format, "IEEE")) { writer = new IEEEIO(filename, IObase::Create); #ifdef HDF4 } else if (CCTK_Equals(out3D_format, "HDF4")) { writer = new HDFIO(filename, IObase::Create); #endif #ifdef HDF5 } else if (CCTK_Equals(out3D_format, "HDF5")) { writer = new H5IO(filename, IObase::Create); #endif } else { assert (0); } delete writer; writer = 0; } } // Open the file if (CCTK_Equals(out3D_format, "IEEE")) { writer = new IEEEIO(filename, IObase::Append); #ifdef HDF4 } else if (CCTK_Equals(out3D_format, "HDF4")) { writer = new HDFIO(filename, IObase::Append); #endif #ifdef HDF5 } else if (CCTK_Equals(out3D_format, "HDF5")) { writer = new H5IO(filename, IObase::Append); #endif } else { assert (0); } assert (writer->isValid()); amrwriter = new AMRwriter(*writer); // Set datatype assert (CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL8 || (sizeof(CCTK_REAL) == sizeof(CCTK_REAL8) && CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL)); // TODO: Set datatype correctly amrwriter->setType (IObase::Float64); const int gpdim = CCTK_GroupDimI(group); // Set coordinate information CCTK_REAL lower[dim], upper[dim]; double origin[dim], delta[dim], timestep; for (int d=0; dcctk_delta_space[d]; } timestep = cgh->cctk_delta_time; amrwriter->setTopLevelParameters (gpdim, origin, delta, timestep, maxreflevels); // Set refinement information int interlevel_timerefinement; int interlevel_spacerefinement[dim]; int initial_gridplacementrefinement[dim]; interlevel_timerefinement = hh->reffact; for (int d=0; dreffact; initial_gridplacementrefinement[d] = 1; } amrwriter->setRefinement (interlevel_timerefinement, interlevel_spacerefinement, initial_gridplacementrefinement); // Set level amrwriter->setLevel (reflevel); // Set current time amrwriter->setTime (cgh->cctk_iteration); } // Traverse all components on this refinement and multigrid level BEGIN_COMPONENT_LOOP(cgh, grouptype) { const ggf* ff = 0; assert (var < (int)arrdata[group].data.size()); ff = (ggf*)arrdata[group].data[var]; CCTK_VInfo (CCTK_THORNSTRING, "bogus reflevel,component,mglevel %d,%d,%d",reflevel,component,mglevel); const gdata* const data = (*ff) (tl, rl, component, mglevel); // Make temporary copy on processor 0 bbox ext = data->extent(); vect lo = ext.lower(); vect hi = ext.upper(); vect str = ext.stride(); // Ignore ghost zones if desired for (int d=0; dcctk_bbox[2*d ] && !out3D_output_outer_boundary) ? -1 : out3D_max_num_lower_ghosts; const int max_upper_ghosts = (cgh->cctk_bbox[2*d+1] && !out3D_output_outer_boundary) ? -1 : out3D_max_num_upper_ghosts; const int num_lower_ghosts = max_lower_ghosts == -1 ? cgh->cctk_nghostzones[d] : min(out3D_max_num_lower_ghosts, cgh->cctk_nghostzones[d]); const int num_upper_ghosts = max_upper_ghosts == -1 ? cgh->cctk_nghostzones[d] : min(out3D_max_num_upper_ghosts, cgh->cctk_nghostzones[d]); lo[d] += (cgh->cctk_nghostzones[d] - num_lower_ghosts) * str[d]; hi[d] -= (cgh->cctk_nghostzones[d] - num_upper_ghosts) * str[d]; } ext = bbox(lo,hi,str); gdata* const tmp = data->make_typed (); tmp->allocate (ext, 0); tmp->copy_from (data, ext); // Write data if (CCTK_MyProc(cgh)==0) { int origin[dim], dims[dim]; for (int d=0; dwrite (origin, dims, (void*)tmp->storage()); } // Delete temporary copy delete tmp; } END_COMPONENT_LOOP; // Close the file if (CCTK_MyProc(cgh)==0) { delete amrwriter; amrwriter = 0; delete writer; writer = 0; } // Don't truncate again do_truncate[n] = false; return 0; } int TimeToOutput (const cGH* const cgh, const int vindex) { DECLARE_CCTK_PARAMETERS; assert (vindex>=0 && vindex<(int)last_output[reflevel].size()); const int myoutevery = GetIntParameter("out3D_every", out_every); if (myoutevery < 0) { // Nothing should be output at all return 0; } if (cgh->cctk_iteration % myoutevery != 0) { // Nothing should be output during this iteration return 0; } if (! CheckForVariable(cgh, GetStringParameter("out3D_vars",""), vindex)) { // This variable should not be output return 0; } if (last_output[reflevel][vindex] == cgh->cctk_iteration) { // Has already been output during this iteration char* varname = CCTK_FullName(vindex); CCTK_VWarn (5, __LINE__, __FILE__, CCTK_THORNSTRING, "Skipping output for variable \"%s\", because this variable " "has already been output during the current iteration -- " "probably via a trigger during the analysis stage", varname); free (varname); return 0; } assert (last_output[reflevel][vindex] < cgh->cctk_iteration); // Should be output during this iteration return 1; } int TriggerOutput (const cGH* const cgh, const int vindex) { assert (vindex>=0 && vindexcctk_iteration; return retval; } int InputGH (const cGH* const cgh) { int retval = 0; for (int vindex=0; vindex=0 && n=0 && group<(int)Carpet::arrdata.size()); const int n0 = CCTK_FirstVarIndexI(group); assert (n0>=0 && n0=0 && varisValid()) { CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING, "Could not open file \"%s\" for reading", filename); } assert (reader->isValid()); if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Reading AMR info"); amrreader = new AmrGridReader(*reader); // Read information about dataset if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Reading dataset info"); IObase::DataType numbertype; reader->readInfo (numbertype, rank, dims); nbytes = IObase::nBytes(numbertype,rank,dims); if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "type=%d rank=%d dims=[%d,%d,%d] nbytes=%d", (int)numbertype, rank, dims[0], dims[1], dims[2], nbytes); // Check rank assert (rank==gpdim); // Check datatype // TODO: Check datatype correctly assert (CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL8 || (sizeof(CCTK_REAL) == sizeof(CCTK_REAL8) && CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL)); // TODO: check grid spacing // Number of datasets if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Reading number of datasets"); ndatasets = reader->nDatasets(); if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "ndatasets=%d", ndatasets); assert (ndatasets>=0); } // Broadcast rank, dimensions, and nbytes MPI_Bcast (&rank, 1, MPI_INT, 0, dist::comm); assert (rank>=1); MPI_Bcast (&dims, rank, MPI_INT, 0, dist::comm); for (int d=0; d=0); MPI_Bcast (&nbytes, 1, MPI_INT, 0, dist::comm); assert (nbytes>=0); // Broadcast number of datasets MPI_Bcast (&ndatasets, 1, MPI_INT, 0, dist::comm); assert (ndatasets>=0); // Read all datasets // TODO: read only some datasets for (int dataset=0; datasetgetGrid(dataset); assert (amrgrid!=0); assert (amrgrid->data!=0); // If iorigin attribute is absent, assume file has unigrid // data. Initialize iorigin to 0. IObase::DataType atype; int alength; if (reader->readAttributeInfo("iorigin", atype, alength) < 0) { for (int d=0; diorigin[d] = 0; } } for (int d=0; diorigin[d]; amr_dims[d] = amrgrid->dims[d]; } for (int d=gpdim; d* ff = 0; assert (var < (int)arrdata[group].data.size()); ff = (ggf*)arrdata[group].data[var]; gdata* const data = (*ff) (tl, rl, component, mglevel); // Create temporary data storage on processor 0 const vect str = vect(reflevelfact); const vect lb = vect(amr_origin) * str; const vect ub = lb + (vect(amr_dims) - 1) * str; const bbox ext(lb,ub,str); gdata* const tmp = data->make_typed (); if (CCTK_MyProc(cgh)==0) { tmp->allocate (ext, 0, amrgrid->data); } else { tmp->allocate (ext, 0, 0); } // Copy into grid function data->copy_from (tmp, ext); // Delete temporary copy delete tmp; } END_COMPONENT_LOOP; if (CCTK_MyProc(cgh)==0) { free (amrgrid->data); free (amrgrid); amrgrid = 0; } } // loop over datasets // Close the file if (CCTK_MyProc(cgh)==0) { if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Deleting AMR info"); delete amrreader; amrreader = 0; if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Closing file"); delete reader; reader = 0; } return 0; } int CarpetIOFlexIOReadData (CCTK_ARGUMENTS) { DECLARE_CCTK_ARGUMENTS; return InputGH(cctkGH); } const char* GetStringParameter (const char* const parametername, const char* const fallback) { if (CCTK_ParameterQueryTimesSet (parametername, CCTK_THORNSTRING) > 0) { int ptype; const char* const* const ppval = (const char* const*)CCTK_ParameterGet (parametername, CCTK_THORNSTRING, &ptype); assert (ppval); const char* const pval = *ppval; assert (ptype == PARAMETER_STRING); return pval; } return fallback; } int GetIntParameter (const char* const parametername, int fallback) { if (CCTK_ParameterQueryTimesSet (parametername, CCTK_THORNSTRING) > 0) { int ptype; const int* const ppval = (const int*)CCTK_ParameterGet (parametername, CCTK_THORNSTRING, &ptype); assert (ppval); const int pval = *ppval; assert (ptype == PARAMETER_INT); return pval; } return fallback; } bool CheckForVariable (const cGH* const cgh, const char* const varlist, const int vindex) { const int numvars = CCTK_NumVars(); assert (vindex>=0 && vindex flags(numvars); CCTK_TraverseString (varlist, SetFlag, &flags, CCTK_GROUP_OR_VAR); return flags[vindex]; } void SetFlag (int index, const char* optstring, void* arg) { vector& flags = *(vector*)arg; flags[index] = true; } } // namespace CarpetIOFlexIO