Page MenuHomec4science

AMEL.cc
No OneTemporary

File Metadata

Created
Sun, Jun 9, 21:16
/* ./amel/AMEL.cpp
**********************************
author : Guillaume ANCIAUX (guillaume.anciaux@epfl.ch, g.anciaux@laposte.net)
The LibMultiScale is a C++ parallel framework for the multiscale
coupling methods dedicated to material simulations. This framework
provides an API which makes it possible to program coupled simulations
and integration of already existing codes.
This Project is done in a collaboration between
EPFL within ENAC-LSMS (http://lsms.epfl.ch/) and
INRIA Bordeaux, ScAlApplix (http://www.labri.fr/projet/scalapplix/).
This software is governed by the CeCILL-C license under French law and
abiding by the rules of distribution of free software. You can use,
modify and/ or redistribute the software under the terms of the CeCILL-C
license as circulated by CEA, CNRS and INRIA at the following URL
"http://www.cecill.info".
As a counterpart to the access to the source code and rights to copy,
modify and redistribute granted by the license, users are provided only
with a limited warranty and the software's author, the holder of the
economic rights, and the successive licensors have only limited
liability.
In this respect, the user's attention is drawn to the risks associated
with loading, using, modifying and/or developing or reproducing the
software by the user in light of its specific status of free software,
that may mean that it is complicated to manipulate, and that also
therefore means that it is reserved for developers and experienced
professionals having in-depth computer knowledge. Users are therefore
encouraged to load and test the software's suitability as regards their
requirements in conditions enabling the security of their systems and/or
data to be ensured and, more generally, to use and operate it in the
same conditions as regards security.
The fact that you are presently reading this means that you have had
knowledge of the CeCILL-C license and that you accept its terms.
***********************************/
/* -------------------------------------------------------------------------- */
#define TIMER
/* -------------------------------------------------------------------------- */
#include "action_interface.hh"
#include "algebraic_parser.hh"
#include "domain_multiscale.hh"
#include "factory_multiscale.hh"
#include "filter_interface.hh"
#include "lm_common.hh"
/* -------------------------------------------------------------------------- */
#include <iomanip>
#include <iostream>
#include <mpi.h>
#include <sstream>
#include <stdexcept>
#include <sys/types.h>
#include <unistd.h>
/* -------------------------------------------------------------------------- */
using namespace libmultiscale;
/* -------------------------------------------------------------------------- */
Real lastdump = -1e30;
Real lastStep = 0;
/* -------------------------------------------------------------------------- */
void printState() {
if (lm_my_proc_id != 0)
return;
struct timespec gtime;
clock_gettime(CLOCK_PROCESS_CPUTIME_ID, &gtime);
if (gtime.tv_sec - lastdump < 2)
return;
Real nstep_done = -lastStep + current_step;
std::stringstream info_steps;
info_steps << std::setw(3) << std::setfill(' ')
<< 100 * current_step / nb_step << "% - "
#ifdef USE_COMPACT_STATUS
<< "step"
#else
<< "passing step "
#endif
<< std::setw(5) << std::setfill(' ') << current_step << "/"
<< std::setw(5) << std::setfill(' ') << nb_step;
if (current_step > 0) {
Real step_per_seconds = nstep_done / (gtime.tv_sec - lastdump);
Real remaining_time = 1. * (nb_step - current_step) / step_per_seconds;
UInt remaining_hours = int(remaining_time / 60 / 60);
UInt remaining_minutes = int(remaining_time / 60 - remaining_hours * 60);
UInt remaining_seconds = int(remaining_time - remaining_minutes * 60 -
remaining_hours * 60 * 60);
info_steps << " " << std::setw(5) << std::setfill(' ') << step_per_seconds
#ifdef USE_COMPACT_STATUS
<< "/sec. ETA: "
#else
<< " steps/seconds "
<< "remaining "
#endif
<< std::setw(5) << std::setfill(' ') << remaining_hours
#ifdef USE_COMPACT_STATUS
<< ":"
#else
<< " hours "
#endif
#ifdef USE_COMPACT_STATUS
<< std::setw(2) << std::setfill('0') << remaining_minutes << ":"
#else
<< std::setw(2) << std::setfill(' ') << remaining_minutes
<< " minutes "
#endif
#ifdef USE_COMPACT_STATUS
<< std::setw(2) << std::setfill('0') << remaining_seconds;
#else
<< std::setw(2) << std::setfill(' ') << remaining_seconds
<< " seconds";
#endif
}
std::string clown = info_steps.str();
// std::cerr << "\r" << clown;
std::cerr << info_steps.str() << std::endl;
lastStep = current_step;
lastdump = gtime.tv_sec;
}
/* -------------------------------------------------------------------------- */
// void echange_args(const string &* argv,char **my_args);
extern std::string lm_release_info;
/* -------------------------------------------------------------------------- */
static void Usage(const std::string &name) {
std::cout
<< "Usage : " << name << " config_file nb_time_step\n"
<< "\t config_file : global config file of the simulation\n"
<< "\t nb_time_step : number of time step wanted to be done (>=0)\n";
std::cout << "*****************************************" << std::endl;
std::cout << "release info" << std::endl;
std::cout << "*****************************************" << std::endl;
std::cout << lm_release_info << std::endl;
lm_exit(LM_EXIT_FAILURE);
}
/* -------------------------------------------------------------------------- */
int main(int argc, char **argv) {
loadModules(argc, argv); // loading the modules macro
if (argc != 3) {
if (lm_my_proc_id == 0)
Usage(argv[0]);
lm_exit(LM_EXIT_FAILURE);
}
nb_step = atoi(argv[2]);
nb_step_next_event = nb_step;
if (nb_step == UINT_MAX)
Usage(argv[0]);
MPI_Barrier(MPI_COMM_WORLD);
DOWAIT_AT_STARTUP;
STARTTIMER("Init0");
DomainMultiScale &dom = DomainMultiScale::getManager();
dom.build(argv[1]);
auto &actions = ActionManager::getManager();
auto &filters = FilterManager::getManager();
STOPTIMER("Init0");
nb_step += current_step;
MPI_Barrier(MPI_COMM_WORLD);
bool shouldPrintState = true;
char *varenv = getenv("PRINT_STATE");
if (varenv)
shouldPrintState = true;
STARTTIMER("Main_loop");
Real time_step = 1e300;
dom.for_each([&time_step](auto &d) {
Real ts = d.template getParam<Quantity<Time>>("TIMESTEP");
time_step = std::min(time_step, ts);
});
DUMP("Multiscale timestep:" << time_step, DBG_MESSAGE);
// Real min_dt = dom.getTimeStep();
// current_time = current_step * min_dt;
#ifdef LIBMULTISCALE_TIMER
UInt freq = 1;
if (nb_step >= 10000)
freq = 100;
#endif
STARTTIMER("updateAcceleration");
dom.for_each([](auto &d) { d.updateGradient(); });
dom.coupling(COUPLING_STEP2);
dom.for_each([](auto &d) { d.updateAcceleration(); });
STOPTIMER("updateAcceleration");
for (; current_step < nb_step; ++current_step, current_time += time_step) {
if (shouldPrintState)
printState();
try {
STARTTIMER("preDump");
current_stage = PRE_DUMP;
actions.action();
STOPTIMER("preDump");
current_stage = PRE_STEP1;
STARTTIMER("preAction1");
actions.action();
STOPTIMER("preAction1");
STARTTIMER("newmarkPredictor");
dom.for_each([time_step](auto &d) {
auto &&v = d.primalTimeDerivative();
auto &&p = d.primal();
auto &&acc = d.acceleration();
v += 0.5 * time_step * acc;
p += time_step * v;
d.changeRelease();
});
STOPTIMER("newmarkpredictor");
STARTTIMER("enforceCompatibility");
dom.for_each([](auto &d) { d.enforceCompatibility(); });
STOPTIMER("enforceCompatibility");
STARTTIMER("CouplingStep1");
dom.coupling(COUPLING_STEP1);
STOPTIMER("CouplingStep1");
current_stage = PRE_STEP2;
STARTTIMER("preAction2");
actions.action();
STOPTIMER("preAction2");
STARTTIMER("updateForce");
dom.for_each([](auto &d) { d.updateGradient(); });
STOPTIMER("updateForce");
STARTTIMER("CouplingStep2");
dom.coupling(COUPLING_STEP2);
STOPTIMER("CouplingStep2");
current_stage = PRE_STEP3;
STARTTIMER("preAction3");
actions.action();
STOPTIMER("preAction3");
STARTTIMER("updateAcceleration");
dom.for_each([](auto &d) { d.updateAcceleration(); });
STOPTIMER("updateAcceleration");
STARTTIMER("newmarkCorrector");
dom.for_each([time_step](auto &d) {
auto &&v = d.primalTimeDerivative();
auto &&acc = d.acceleration();
v += 0.5 * time_step * acc;
d.changeRelease();
});
STOPTIMER("newmarkCorrector");
STARTTIMER("CouplingStep3");
dom.coupling(COUPLING_STEP3);
STOPTIMER("CouplingStep3");
current_stage = PRE_STEP4;
STARTTIMER("preAction4");
actions.action();
STOPTIMER("preAction4");
STARTTIMER("CouplingStep4");
dom.coupling(COUPLING_STEP4);
STOPTIMER("CouplingStep4");
#ifdef LIBMULTISCALE_TIMER
if (current_step % freq == 0)
dumpTimesInColumn(lm_my_proc_id, current_step);
#endif
} catch (LibMultiScaleException &e) {
LM_FATAL_RE(e, " at stage " << current_stage << std::endl);
} catch (std::runtime_error &err) {
std::stringstream error_stream;
error_stream << " at stage " << current_stage << ": Caught '"
<< err.what() << "'";
LM_FATAL(error_stream.str());
}
}
current_stage = PRE_DUMP;
actions.action();
current_stage = PRE_STEP1;
actions.action();
STOPTIMER("Main_loop");
dom.destroy();
actions.destroy();
filters.destroy();
closeModules(); // closing the modules macro
std::cerr << std::endl;
lm_exit(LM_EXIT_SUCCESS);
}
/* --------------------------------------------------------------------------
*/

Event Timeline