| 1 | /*
|
|---|
| 2 | A simple 2D hydro code
|
|---|
| 3 | (C) Romain Teyssier : CEA/IRFU -- original F90 code
|
|---|
| 4 | (C) Pierre-Francois Lavallee : IDRIS -- original F90 code
|
|---|
| 5 | (C) Guillaume Colin de Verdiere : CEA/DAM -- for the C version
|
|---|
| 6 | */
|
|---|
| 7 | /*
|
|---|
| 8 |
|
|---|
| 9 | This software is governed by the CeCILL license under French law and
|
|---|
| 10 | abiding by the rules of distribution of free software. You can use,
|
|---|
| 11 | modify and/ or redistribute the software under the terms of the CeCILL
|
|---|
| 12 | license as circulated by CEA, CNRS and INRIA at the following URL
|
|---|
| 13 | "http://www.cecill.info".
|
|---|
| 14 |
|
|---|
| 15 | As a counterpart to the access to the source code and rights to copy,
|
|---|
| 16 | modify and redistribute granted by the license, users are provided only
|
|---|
| 17 | with a limited warranty and the software's author, the holder of the
|
|---|
| 18 | economic rights, and the successive licensors have only limited
|
|---|
| 19 | liability.
|
|---|
| 20 |
|
|---|
| 21 | In this respect, the user's attention is drawn to the risks associated
|
|---|
| 22 | with loading, using, modifying and/or developing or reproducing the
|
|---|
| 23 | software by the user in light of its specific status of free software,
|
|---|
| 24 | that may mean that it is complicated to manipulate, and that also
|
|---|
| 25 | therefore means that it is reserved for developers and experienced
|
|---|
| 26 | professionals having in-depth computer knowledge. Users are therefore
|
|---|
| 27 | encouraged to load and test the software's suitability as regards their
|
|---|
| 28 | requirements in conditions enabling the security of their systems and/or
|
|---|
| 29 | data to be ensured and, more generally, to use and operate it in the
|
|---|
| 30 | same conditions as regards security.
|
|---|
| 31 |
|
|---|
| 32 | The fact that you are presently reading this means that you have had
|
|---|
| 33 | knowledge of the CeCILL license and that you accept its terms.
|
|---|
| 34 |
|
|---|
| 35 | */
|
|---|
| 36 |
|
|---|
| 37 | #include <stdlib.h>
|
|---|
| 38 | #include <unistd.h>
|
|---|
| 39 | #include <math.h>
|
|---|
| 40 | #include <stdio.h>
|
|---|
| 41 | #include <strings.h>
|
|---|
| 42 | #include <string.h>
|
|---|
| 43 |
|
|---|
| 44 | #include "hmpp.h"
|
|---|
| 45 | #include "parametres.h"
|
|---|
| 46 | #include "hydro_godunov.h"
|
|---|
| 47 | #include "hydro_funcs.h"
|
|---|
| 48 | #include "cclock.h"
|
|---|
| 49 | #include "utils.h"
|
|---|
| 50 | #include "make_boundary.h"
|
|---|
| 51 |
|
|---|
| 52 | #include "riemann.h"
|
|---|
| 53 | #include "qleftright.h"
|
|---|
| 54 | #include "trace.h"
|
|---|
| 55 | #include "slope.h"
|
|---|
| 56 | #include "equation_of_state.h"
|
|---|
| 57 | #include "constoprim.h"
|
|---|
| 58 |
|
|---|
| 59 | #include "cmpflx.h"
|
|---|
| 60 | #include "conservar.h"
|
|---|
| 61 |
|
|---|
| 62 | // variables auxiliaires pour mettre en place le mode resident de HMPP
|
|---|
| 63 | void
|
|---|
| 64 | hydro_godunov(int idimStart, real_t dt, const hydroparam_t H, hydrovar_t * Hv, hydrowork_t * Hw, hydrovarwork_t * Hvw) {
|
|---|
| 65 | // Local variables
|
|---|
| 66 | struct timespec start, end;
|
|---|
| 67 | int j;
|
|---|
| 68 | real_t dtdx;
|
|---|
| 69 | int clear=0;
|
|---|
| 70 |
|
|---|
| 71 | real_t (*e)[H.nxyt];
|
|---|
| 72 | real_t (*flux)[H.nxystep][H.nxyt];
|
|---|
| 73 | real_t (*qleft)[H.nxystep][H.nxyt];
|
|---|
| 74 | real_t (*qright)[H.nxystep][H.nxyt];
|
|---|
| 75 | real_t (*c)[H.nxyt];
|
|---|
| 76 | real_t *uold;
|
|---|
| 77 | int (*sgnm)[H.nxyt];
|
|---|
| 78 | real_t (*qgdnv)[H.nxystep][H.nxyt];
|
|---|
| 79 | real_t (*u)[H.nxystep][H.nxyt];
|
|---|
| 80 | real_t (*qxm)[H.nxystep][H.nxyt];
|
|---|
| 81 | real_t (*qxp)[H.nxystep][H.nxyt];
|
|---|
| 82 | real_t (*q)[H.nxystep][H.nxyt];
|
|---|
| 83 | real_t (*dq)[H.nxystep][H.nxyt];
|
|---|
| 84 |
|
|---|
| 85 | static FILE *fic = NULL;
|
|---|
| 86 |
|
|---|
| 87 | if (fic == NULL && H.prt == 1) {
|
|---|
| 88 | char logname[256];
|
|---|
| 89 | sprintf(logname, "TRACE.%04d_%04d.txt", H.nproc, H.mype);
|
|---|
| 90 | fic = fopen(logname, "w");
|
|---|
| 91 | }
|
|---|
| 92 |
|
|---|
| 93 | WHERE("hydro_godunov");
|
|---|
| 94 |
|
|---|
| 95 | // int hmppGuard = 1;
|
|---|
| 96 | int idimIndex = 0;
|
|---|
| 97 |
|
|---|
| 98 | for (idimIndex = 0; idimIndex < 2; idimIndex++) {
|
|---|
| 99 | int idim = (idimStart - 1 + idimIndex) % 2 + 1;
|
|---|
| 100 | // constant
|
|---|
| 101 | dtdx = dt / H.dx;
|
|---|
| 102 |
|
|---|
| 103 | // Update boundary conditions
|
|---|
| 104 | if (H.prt) {
|
|---|
| 105 | fprintf(fic, "godunov %d %le %le\n", idim, dt, H.t);
|
|---|
| 106 | PRINTUOLD(fic, H, Hv);
|
|---|
| 107 | }
|
|---|
| 108 | // if (H.mype == 1) fprintf(fic, "Hydro makes boundary.\n");
|
|---|
| 109 | start = cclock();
|
|---|
| 110 | make_boundary(idim, H, Hv);
|
|---|
| 111 | end = cclock();
|
|---|
| 112 | functim[TIM_MAKBOU] += ccelaps(start, end);
|
|---|
| 113 |
|
|---|
| 114 | if (H.prt) {fprintf(fic, "MakeBoundary\n");}
|
|---|
| 115 | PRINTUOLD(fic, H, Hv);
|
|---|
| 116 |
|
|---|
| 117 | uold = Hv->uold;
|
|---|
| 118 | qgdnv = (real_t (*)[H.nxystep][H.nxyt]) Hvw->qgdnv;
|
|---|
| 119 | flux = (real_t (*)[H.nxystep][H.nxyt]) Hvw->flux;
|
|---|
| 120 | c = (real_t (*)[H.nxyt]) Hw->c;
|
|---|
| 121 | e = (real_t (*)[H.nxyt]) Hw->e;
|
|---|
| 122 | qleft = (real_t (*)[H.nxystep][H.nxyt]) Hvw->qleft;
|
|---|
| 123 | qright = (real_t (*)[H.nxystep][H.nxyt]) Hvw->qright;
|
|---|
| 124 | sgnm = (int (*)[H.nxyt]) Hw->sgnm;
|
|---|
| 125 | q = (real_t (*)[H.nxystep][H.nxyt]) Hvw->q;
|
|---|
| 126 | dq = (real_t (*)[H.nxystep][H.nxyt]) Hvw->dq;
|
|---|
| 127 | u = (real_t (*)[H.nxystep][H.nxyt]) Hvw->u;
|
|---|
| 128 | qxm = (real_t (*)[H.nxystep][H.nxyt]) Hvw->qxm;
|
|---|
| 129 | qxp = (real_t (*)[H.nxystep][H.nxyt]) Hvw->qxp;
|
|---|
| 130 |
|
|---|
| 131 | int Hmin, Hmax, Hstep;
|
|---|
| 132 | int Hdimsize;
|
|---|
| 133 | int Hndim_1;
|
|---|
| 134 |
|
|---|
| 135 | if (idim == 1) {
|
|---|
| 136 | Hmin = H.jmin + ExtraLayer;
|
|---|
| 137 | Hmax = H.jmax - ExtraLayer;
|
|---|
| 138 | Hdimsize = H.nxt;
|
|---|
| 139 | Hndim_1 = H.nx + 1;
|
|---|
| 140 | Hstep = H.nxystep;
|
|---|
| 141 | } else {
|
|---|
| 142 | Hmin = H.imin + ExtraLayer;
|
|---|
| 143 | Hmax = H.imax - ExtraLayer;
|
|---|
| 144 | Hdimsize = H.nyt;
|
|---|
| 145 | Hndim_1 = H.ny + 1;
|
|---|
| 146 | Hstep = H.nxystep;
|
|---|
| 147 | }
|
|---|
| 148 |
|
|---|
| 149 | if (!H.nstep && idim == 1) {
|
|---|
| 150 | /* LM -- HERE a more secure implementation should be used: a new parameter ? */
|
|---|
| 151 | }
|
|---|
| 152 | // if (H.mype == 1) fprintf(fic, "Hydro computes slices.\n");
|
|---|
| 153 | for (j = Hmin; j < Hmax; j += Hstep) {
|
|---|
| 154 | // we try to compute many slices each pass
|
|---|
| 155 | int jend = j + Hstep;
|
|---|
| 156 | if (jend >= Hmax)
|
|---|
| 157 | jend = Hmax;
|
|---|
| 158 | int slices = jend - j; // numbre of slices to compute
|
|---|
| 159 | // fprintf(stderr, "Godunov idim=%d, j=%d %d \n", idim, j, slices);
|
|---|
| 160 |
|
|---|
| 161 | if (clear) Dmemset((H.nxyt) * H.nxystep * H.nvar, (real_t *) dq, 0);
|
|---|
| 162 | start = cclock();
|
|---|
| 163 | gatherConservativeVars(idim, j, H.imin, H.imax, H.jmin, H.jmax, H.nvar, H.nxt, H.nyt, H.nxyt, slices, Hstep, uold,
|
|---|
| 164 | u);
|
|---|
| 165 | end = cclock();
|
|---|
| 166 | functim[TIM_GATCON] += ccelaps(start, end);
|
|---|
| 167 | if (H.prt) {fprintf(fic, "ConservativeVars %d %d %d %d %d %d\n", H.nvar, H.nxt, H.nyt, H.nxyt, slices, Hstep);}
|
|---|
| 168 | PRINTARRAYV2(fic, u, Hdimsize, "u", H);
|
|---|
| 169 |
|
|---|
| 170 | if (clear) Dmemset((H.nxyt) * H.nxystep * H.nvar, (real_t *) dq, 0);
|
|---|
| 171 |
|
|---|
| 172 | // Convert to primitive variables
|
|---|
| 173 | start = cclock();
|
|---|
| 174 | constoprim(Hdimsize, H.nxyt, H.nvar, H.smallr, slices, Hstep, u, q, e);
|
|---|
| 175 | end = cclock();
|
|---|
| 176 | functim[TIM_CONPRI] += ccelaps(start, end);
|
|---|
| 177 | PRINTARRAY(fic, e, Hdimsize, "e", H);
|
|---|
| 178 | PRINTARRAYV2(fic, q, Hdimsize, "q", H);
|
|---|
| 179 |
|
|---|
| 180 | start = cclock();
|
|---|
| 181 | equation_of_state(0, Hdimsize, H.nxyt, H.nvar, H.smallc, H.gamma, slices, Hstep, e, q, c);
|
|---|
| 182 | end = cclock();
|
|---|
| 183 | functim[TIM_EOS] += ccelaps(start, end);
|
|---|
| 184 | PRINTARRAY(fic, c, Hdimsize, "c", H);
|
|---|
| 185 | PRINTARRAYV2(fic, q, Hdimsize, "q", H);
|
|---|
| 186 |
|
|---|
| 187 | // Characteristic tracing
|
|---|
| 188 | if (H.iorder != 1) {
|
|---|
| 189 | start = cclock();
|
|---|
| 190 | slope(Hdimsize, H.nvar, H.nxyt, H.slope_type, slices, Hstep, q, dq);
|
|---|
| 191 | end = cclock();
|
|---|
| 192 | functim[TIM_SLOPE] += ccelaps(start, end);
|
|---|
| 193 | PRINTARRAYV2(fic, dq, Hdimsize, "dq", H);
|
|---|
| 194 | }
|
|---|
| 195 |
|
|---|
| 196 | if (clear) Dmemset(H.nxyt * H.nxystep * H.nvar, (real_t *) qxm, 0);
|
|---|
| 197 | if (clear) Dmemset(H.nxyt * H.nxystep * H.nvar, (real_t *) qxp, 0);
|
|---|
| 198 | if (clear) Dmemset(H.nxyt * H.nxystep * H.nvar, (real_t *) qleft, 0);
|
|---|
| 199 | if (clear) Dmemset(H.nxyt * H.nxystep * H.nvar, (real_t *) qright, 0);
|
|---|
| 200 | if (clear) Dmemset(H.nxyt * H.nxystep * H.nvar, (real_t *) flux, 0);
|
|---|
| 201 | if (clear) Dmemset(H.nxyt * H.nxystep * H.nvar, (real_t *) qgdnv, 0);
|
|---|
| 202 | start = cclock();
|
|---|
| 203 | trace(dtdx, Hdimsize, H.scheme, H.nvar, H.nxyt, slices, Hstep, q, dq, c, qxm, qxp);
|
|---|
| 204 | end = cclock();
|
|---|
| 205 | functim[TIM_TRACE] += ccelaps(start, end);
|
|---|
| 206 | PRINTARRAYV2(fic, qxm, Hdimsize, "qxm", H);
|
|---|
| 207 | PRINTARRAYV2(fic, qxp, Hdimsize, "qxp", H);
|
|---|
| 208 |
|
|---|
| 209 | start = cclock();
|
|---|
| 210 | qleftright(idim, H.nx, H.ny, H.nxyt, H.nvar, slices, Hstep, qxm, qxp, qleft, qright);
|
|---|
| 211 | end = cclock();
|
|---|
| 212 | functim[TIM_QLEFTR] += ccelaps(start, end);
|
|---|
| 213 | PRINTARRAYV2(fic, qleft, Hdimsize, "qleft", H);
|
|---|
| 214 | PRINTARRAYV2(fic, qright, Hdimsize, "qright", H);
|
|---|
| 215 |
|
|---|
| 216 | start = cclock();
|
|---|
| 217 | riemann(Hndim_1, H.smallr, H.smallc, H.gamma, H.niter_riemann, H.nvar, H.nxyt, slices, Hstep, qleft, qright, qgdnv, sgnm, Hw);
|
|---|
| 218 | end = cclock();
|
|---|
| 219 | functim[TIM_RIEMAN] += ccelaps(start, end);
|
|---|
| 220 | PRINTARRAYV2(fic, qgdnv, Hdimsize, "qgdnv", H);
|
|---|
| 221 |
|
|---|
| 222 | start = cclock();
|
|---|
| 223 | cmpflx(Hdimsize, H.nxyt, H.nvar, H.gamma, slices, Hstep, qgdnv, flux);
|
|---|
| 224 | end = cclock();
|
|---|
| 225 | functim[TIM_CMPFLX] += ccelaps(start, end);
|
|---|
| 226 | PRINTARRAYV2(fic, flux, Hdimsize, "flux", H);
|
|---|
| 227 | PRINTARRAYV2(fic, u, Hdimsize, "u", H);
|
|---|
| 228 |
|
|---|
| 229 | start = cclock();
|
|---|
| 230 | updateConservativeVars(idim, j, dtdx, H.imin, H.imax, H.jmin, H.jmax, H.nvar, H.nxt, H.nyt, H.nxyt, slices, Hstep,
|
|---|
| 231 | uold, u, flux);
|
|---|
| 232 | end = cclock();
|
|---|
| 233 | functim[TIM_UPDCON] += ccelaps(start, end);
|
|---|
| 234 | PRINTUOLD(fic, H, Hv);
|
|---|
| 235 | } // for j
|
|---|
| 236 |
|
|---|
| 237 | if (H.prt) {
|
|---|
| 238 | fprintf(fic, "[%d] After pass %d\n", H.mype, idim);
|
|---|
| 239 | PRINTUOLD(fic, H, Hv);
|
|---|
| 240 | }
|
|---|
| 241 | }
|
|---|
| 242 |
|
|---|
| 243 | if ((H.t + dt >= H.tend) || (H.nstep + 1 >= H.nstepmax)) {
|
|---|
| 244 | /* LM -- HERE a more secure implementation should be used: a new parameter ? */
|
|---|
| 245 | }
|
|---|
| 246 |
|
|---|
| 247 | } // hydro_godunov
|
|---|
| 248 |
|
|---|
| 249 |
|
|---|
| 250 | // EOF
|
|---|