Esempio n. 1
0
/** Replaces observation errors in the observation file with the modified
 * values. The original values are stored as "std_orig".
 */
void das_addmodifiederrors(dasystem* das, char fname[])
{
    int ncid;
    int dimid_nobs[1];
    size_t nobs;
    int varid_std;
    double* std;
    int i;

    if (rank != 0)
        return;

    ncw_open(fname, NC_WRITE, &ncid);
    ncw_inq_dimid(fname, ncid, "nobs", dimid_nobs);
    ncw_inq_dimlen(fname, ncid, dimid_nobs[0], &nobs);

    ncw_redef(fname, ncid);
    if (ncw_var_exists(ncid, "std_orig"))
        enkf_quit("\"observations.nc\" has already been modified by `enkf_calc\'. To proceed please remove observations*.nc and rerun `enkf_prep\'.");
    ncw_rename_var(fname, ncid, "std", "std_orig");
    ncw_def_var(fname, ncid, "std", NC_FLOAT, 1, dimid_nobs, &varid_std);
    ncw_enddef(fname, ncid);

    std = malloc(nobs * sizeof(double));

    for (i = 0; i < (int) nobs; ++i)
        std[i] = das->obs->data[i].std;
    ncw_put_var_double(fname, ncid, varid_std, std);

    free(std);

    ncw_close(fname, ncid);
}
Esempio n. 2
0
/** Adds forecast observations and forecast ensemble spread to the observation
 * file.
 */
void das_addforecast(dasystem* das, char fname[])
{
    int ncid;
    int dimid_nobs[1];
    size_t nobs;
    int varid_Hx, varid_spread;
    double* Hx;
    int o;

    if (das->obs->nobs == 0)
        return;
    if (rank != 0)
        return;

    assert(das->s_mode == S_MODE_HA_f);

    ncw_open(fname, NC_WRITE, &ncid);
    if (ncw_var_exists(ncid, "Hx_f")) {
        enkf_printf("  Hx already added to \"%s\" (skipping)\n", fname);
        goto finish;
    }

    ncw_inq_dimid(fname, ncid, "nobs", dimid_nobs);
    ncw_inq_dimlen(fname, ncid, dimid_nobs[0], &nobs);
    ncw_redef(fname, ncid);
    ncw_def_var(fname, ncid, "Hx_f", NC_FLOAT, 1, dimid_nobs, &varid_Hx);
    ncw_def_var(fname, ncid, "std_f", NC_FLOAT, 1, dimid_nobs, &varid_spread);
    ncw_enddef(fname, ncid);

    Hx = calloc(nobs, sizeof(double));
    for (o = 0; o < (int) nobs; ++o)
        Hx[o] = das->obs->data[o].value - das->s_f[o];

    ncw_put_var_double(fname, ncid, varid_Hx, Hx);
    ncw_put_var_double(fname, ncid, varid_spread, das->std_f);

    free(Hx);

  finish:
    ncw_close(fname, ncid);
}
Esempio n. 3
0
/** Add analysed observation estimates and ensemble spread to FNAME_SOBS.
 */
void das_addanalysis(dasystem* das, char fname[])
{
    int ncid;
    int dimid_nobs[1];
    size_t nobs;
    int varid_Hx, varid_spread;
    observation* data = das->obs->data;
    double* s;
    int i;

    if (rank != 0)
        return;

    assert(das->s_mode == S_MODE_HA_a);

    ncw_open(fname, NC_WRITE, &ncid);
    ncw_inq_dimid(fname, ncid, "nobs", dimid_nobs);
    ncw_inq_dimlen(fname, ncid, dimid_nobs[0], &nobs);

    ncw_redef(fname, ncid);
    ncw_def_var(fname, ncid, "Hx_a", NC_FLOAT, 1, dimid_nobs, &varid_Hx);
    ncw_def_var(fname, ncid, "std_a", NC_FLOAT, 1, dimid_nobs, &varid_spread);
    ncw_enddef(fname, ncid);

    ncw_put_var_double(fname, ncid, varid_spread, das->std_a);
    s = malloc(nobs * sizeof(double));
    /*
     * the obs are still sorted by ij 
     */
    for (i = 0; i < (int) nobs; ++i)
        s[data[i].id] = data[i].value;
    for (i = 0; i < (int) nobs; ++i)
        s[i] -= das->s_a[i];
    ncw_put_var_double(fname, ncid, varid_Hx, s);
    free(s);

    ncw_close(fname, ncid);
}
Esempio n. 4
0
/** Replaces observation errors in the observation file with the modified
 * values. The original values are stored as "std_orig".
 */
void das_addmodifiederrors(dasystem* das, char fname[])
{
    int ncid;
    int dimid_nobs[1];
    size_t nobs;
    int varid_std;
    double* std;
    int i;
    double da_julday = NaN;

    if (rank != 0)
        return;

    ncw_open(fname, NC_WRITE, &ncid);
    ncw_inq_dimid(ncid, "nobs", dimid_nobs);
    ncw_inq_dimlen(ncid, dimid_nobs[0], &nobs);

    ncw_get_att_double(ncid, NC_GLOBAL, "DA_JULDAY", &da_julday);
    if (!enkf_noobsdatecheck && (isnan(da_julday) || fabs(das->obs->da_date - da_julday) > 1e-6))
        enkf_quit("\"observations.nc\" from a different cycle");

    if (ncw_var_exists(ncid, "std_orig")) {
        enkf_printf("    nothing to do\n");
        ncw_close(ncid);
        return;
    }

    ncw_redef(ncid);
    ncw_rename_var(ncid, "std", "std_orig");
    ncw_def_var(ncid, "std", NC_FLOAT, 1, dimid_nobs, &varid_std);
    ncw_enddef(ncid);

    std = malloc(nobs * sizeof(double));

    for (i = 0; i < (int) nobs; ++i)
        std[i] = das->obs->data[i].std;
    ncw_put_var_double(ncid, varid_std, std);

    free(std);

    ncw_close(ncid);
}
Esempio n. 5
0
void das_getHE(dasystem* das)
{
    observations* obs = das->obs;
    model* m = das->m;
    ENSOBSTYPE* Hx = NULL;
    int i, e;

    das->s_mode = S_MODE_HE_f;
    if (obs->nobs == 0)
        return;

    if (das->nmem <= 0)
        das_getnmem(das);
    enkf_printf("    ensemble size = %d\n", das->nmem);
    assert(das->nmem > 0);

    distribute_iterations(0, das->nmem - 1, nprocesses, rank, "    ");

    /*
     * ensemble observation array to be filled 
     */
    assert(das->S == NULL);
    das->S = alloc2d(das->nmem, obs->nobs, sizeof(ENSOBSTYPE));
    if (das->mode == MODE_ENOI)
        Hx = calloc(obs->nobs, sizeof(ENSOBSTYPE));

    for (i = 0; i < obs->nobstypes; ++i) {
        obstype* ot = &obs->obstypes[i];
        float*** vvv = NULL;
        float** vv = NULL;
        H_fn H = NULL;
        int mvid;
        int ni, nj, nk;
        int nobs;
        int* obsids;
        char fname[MAXSTRLEN];

        enkf_printf("    %s ", ot->name);
        fflush(stdout);

        mvid = model_getvarid(m, obs->obstypes[i].varname);
        if (mvid < 0)
            enkf_quit("variable \"%s\" required for observation type \"%s\" is not defined", obs->obstypes[i].varname, ot->name);
        if (ot->issurface) {
            model_getvardims(m, mvid, &ni, &nj, NULL);
            vv = alloc2d(nj, ni, sizeof(float));
        } else {
            model_getvardims(m, mvid, &ni, &nj, &nk);
            vvv = alloc3d(nk, nj, ni, sizeof(float));
        }

        /*
         * set H
         */
        H = getH(ot->name, ot->hfunction);

        if (ot->isasync) {
            int t1 = get_tshift(ot->date_min, ot->async_tstep);
            int t2 = get_tshift(ot->date_max, ot->async_tstep);
            int t;

            for (t = t1; t <= t2; ++t) {
                enkf_printf("|");
                obs_find_bytypeandtime(obs, i, t, &nobs, &obsids);
                if (nobs == 0)
                    continue;

                if (das->mode == MODE_ENKF || !enkf_fstatsonly) {
                    for (e = my_first_iteration; e <= my_last_iteration; ++e) {
                        int success = model_getmemberfname_async(m, das->ensdir, ot->varname, ot->name, e + 1, t, fname);

                        H(das, nobs, obsids, fname, e + 1, t, ot->varname, ot->varname2, (ot->issurface) ? (void*) vv : (void*) vvv, das->S[e]);
                        enkf_printf((success) ? "a" : "s");
                        fflush(stdout);
                    }
                }

                if (das->mode == MODE_ENOI) {
                    if (enkf_obstype == OBSTYPE_VALUE) {
                        int success = model_getbgfname_async(m, das->bgdir, ot->varname, ot->name, t, fname);

                        H(das, nobs, obsids, fname, -1, t, ot->varname, ot->varname2, (ot->issurface) ? (void*) vv : (void*) vvv, Hx);
                        enkf_printf((success) ? "A" : "S");
                        fflush(stdout);
                    } else if (enkf_obstype == OBSTYPE_INNOVATION) {
                        Hx[0] = 0;
                        enkf_printf("-");
                        fflush(stdout);
                    }
                }

                free(obsids);
            }
        } else {
            obs_find_bytype(obs, i, &nobs, &obsids);
            if (nobs == 0)
                goto next;

            if (das->mode == MODE_ENKF || !enkf_fstatsonly) {
                for (e = my_first_iteration; e <= my_last_iteration; ++e) {
                    model_getmemberfname(m, das->ensdir, ot->varname, e + 1, fname);
                    H(das, nobs, obsids, fname, e + 1, MAXINT, ot->varname, ot->varname2, (ot->issurface) ? (void*) vv : (void*) vvv, das->S[e]);
                    enkf_printf(".");
                    fflush(stdout);
                }
            }

            if (das->mode == MODE_ENOI) {
                if (enkf_obstype == OBSTYPE_VALUE) {
                    model_getbgfname(m, das->bgdir, ot->varname, fname);
                    H(das, nobs, obsids, fname, -1, MAXINT, ot->varname, ot->varname2, (ot->issurface) ? (void*) vv : (void*) vvv, Hx);
                    enkf_printf("+");
                    fflush(stdout);
                } else if (enkf_obstype == OBSTYPE_INNOVATION) {
                    Hx[0] = 0;
                    enkf_printf("-");
                    fflush(stdout);
                }
            }

            free(obsids);
        }

      next:

        if (ot->issurface)
            free2d(vv);
        else
            free3d(vvv);
        enkf_printf("\n");
    }                           /* for i (over obstypes) */

#if defined(MPI)
    if (das->mode == MODE_ENKF || !enkf_fstatsonly) {
#if !defined(HE_VIAFILE)
        /*
         * communicate HE via MPI
         */
        int ierror, count;

        /*
         * Blocking communications can create a bottleneck in instances with
         * large number of observations (e.g., 2e6 obs., 144 members, 48
         * processes), but asynchronous send/receive seem to work well
         */
        if (rank > 0) {
            MPI_Request request;

            /*
             * send ensemble observations to master
             */
            count = (my_last_iteration - my_first_iteration + 1) * obs->nobs;
            ierror = MPI_Isend(das->S[my_first_iteration], count, MPIENSOBSTYPE, 0, rank, MPI_COMM_WORLD, &request);
            assert(ierror == MPI_SUCCESS);
        } else {
            int r;
            MPI_Request* requests = malloc((nprocesses - 1) * sizeof(MPI_Request));

            /*
             * collect ensemble observations from slaves
             */
            for (r = 1; r < nprocesses; ++r) {
                count = (last_iteration[r] - first_iteration[r] + 1) * obs->nobs;
                ierror = MPI_Irecv(das->S[first_iteration[r]], count, MPIENSOBSTYPE, r, r, MPI_COMM_WORLD, &requests[r - 1]);
                assert(ierror == MPI_SUCCESS);
            }
            ierror = MPI_Waitall(nprocesses - 1, requests, MPI_STATUS_IGNORE);
            assert(ierror == MPI_SUCCESS);
            free(requests);
        }
        /*
         * now send the full set of ensemble observations to slaves
         */
        count = das->nmem * obs->nobs;
        ierror = MPI_Bcast(das->S[0], count, MPIENSOBSTYPE, 0, MPI_COMM_WORLD);
        assert(ierror == MPI_SUCCESS);
#else
        /*
         * communicate HE via file
         */
        {
            int ncid;
            int varid;
            size_t start[2], count[2];

            if (rank == 0) {
                int dimids[2];

                ncw_create(FNAME_HE, NC_CLOBBER | NC_64BIT_OFFSET, &ncid);
                ncw_def_dim(FNAME_HE, ncid, "m", das->nmem, &dimids[0]);
                ncw_def_dim(FNAME_HE, ncid, "p", obs->nobs, &dimids[1]);
                ncw_def_var(FNAME_HE, ncid, "HE", NC_FLOAT, 2, dimids, &varid);
                ncw_close(FNAME_HE, ncid);
            }
            MPI_Barrier(MPI_COMM_WORLD);

            ncw_open(FNAME_HE, NC_WRITE, &ncid);
            ncw_inq_varid(FNAME_HE, ncid, "HE", &varid);
            start[0] = my_first_iteration;
            start[1] = 0;
            count[0] = my_last_iteration - my_first_iteration + 1;
            count[1] = obs->nobs;
            ncw_put_vara_float(FNAME_HE, ncid, varid, start, count, das->S[my_first_iteration]);
            ncw_close(FNAME_HE, ncid);
            MPI_Barrier(MPI_COMM_WORLD);

            ncw_open(FNAME_HE, NC_NOWRITE, &ncid);
            ncw_inq_varid(FNAME_HE, ncid, "HE", &varid);
            ncw_get_var_float(FNAME_HE, ncid, varid, das->S[0]);
            ncw_close(FNAME_HE, ncid);
        }
#endif
    }
#endif

    if (das->mode == MODE_ENOI) {
        /*
         * subtract ensemble mean; add background
         */
        if (!enkf_fstatsonly) {
            double* ensmean = calloc(obs->nobs, sizeof(double));

            for (e = 0; e < das->nmem; ++e) {
                ENSOBSTYPE* Se = das->S[e];

                for (i = 0; i < obs->nobs; ++i)
                    ensmean[i] += Se[i];
            }
            for (i = 0; i < obs->nobs; ++i)
                ensmean[i] /= (double) das->nmem;

            for (e = 0; e < das->nmem; ++e) {
                ENSOBSTYPE* Se = das->S[e];

                for (i = 0; i < obs->nobs; ++i)
                    Se[i] += Hx[i] - ensmean[i];
            }

            free(ensmean);
        } else {
            for (e = 0; e < das->nmem; ++e) {
                ENSOBSTYPE* Se = das->S[e];

                for (i = 0; i < obs->nobs; ++i)
                    Se[i] = Hx[i];
            }
        }
    }

    if (das->mode == MODE_ENOI)
        free(Hx);
}
Esempio n. 6
0
File: nccat.c Progetto: sakov/ncw-c
int main(int argc, char** argv)
{
    int ndims = 0;
    char** dims = NULL;
    int nsrc = 0;
    char** src = NULL;
    char* dst = NULL;
    int i;

    int** dimids = NULL;
    int nvars = 0;
    int* varids = NULL;
    size_t* size = NULL;
    int nvars0;

    if (argc == 1) {
        usage();
        description();
        exit(0);
    }

    i = 1;
    while (i < argc) {
        if (argv[i][0] != '-') {
            usage();
            exit(1);
        }
        if (argv[i][1] == 'd') {
            i++;
            while (i < argc && argv[i][0] != '-') {
                if (ndims % INC == 0)
                    dims = realloc(dims, (ndims + INC) * sizeof(char*));
                dims[ndims] = argv[i];
                ndims++;
                i++;
            }
        } else if (argv[i][1] == 'i') {
            i++;
            while (i < argc && argv[i][0] != '-') {
                if (nsrc % INC == 0)
                    src = realloc(src, (nsrc + INC) * sizeof(char*));
                src[nsrc] = strdup(argv[i]);
                nsrc++;
                i++;
            }
        } else if (argv[i][1] == 'o') {
            i++;
            if (i < argc && argv[i][0] != '-') {
                dst = argv[i];
                i++;
            }
        } else if (argv[i][1] == 'v') {
            verbose = 1;
            i++;
        } else {
            usage();
            exit(1);
        }
    }
    if (argc == 2 && verbose) {
        printf("  ncw version %s\n", ncw_version);
        exit(1);
    }
    if (ndims == 0 || nsrc < 1 || dst == NULL) {
        usage();
        exit(1);
    }

    printf_v("  src:\n");
    for (i = 0; i < nsrc; ++i)
        printf_v("    \"%s\"\n", src[i]);
    printf_v("  dst = \"%s\"\n", dst);
    printf_v("  merged dimensions:\n");
    for (i = 0; i < ndims; ++i)
        printf_v("    %s\n", dims[i]);

    dimids = alloc2d(nsrc, ndims, sizeof(int));
    for (i = 0; i < nsrc; ++i) {
        int ncid, j;

        ncw_open(src[i], NC_NOWRITE, &ncid);
        for (j = 0; j < ndims; ++j)
            ncw_inq_dimid(ncid, dims[j], &dimids[i][j]);
        ncw_close(ncid);
    }

    {
        int ncid;

        ncw_open(src[0], NC_NOWRITE, &ncid);
        ncw_inq_nvars(ncid, &nvars0);
        varids = malloc(nvars0 * sizeof(int));
        printf_v("  merged variables:\n");
        for (i = 0; i < nvars0; ++i) {
            int dimids_now[NC_MAX_VAR_DIMS];
            char varname[NC_MAX_NAME] = "";
            int j;

            ncw_inq_vardimid(ncid, i, dimids_now);
            for (j = 0; j < ndims; ++j) {
                if (dimids_now[0] == dimids[0][j]) {
                    varids[nvars] = i;
                    nvars++;
                    ncw_inq_varname(ncid, i, varname);
                    printf_v("    %s\n", varname);
                    break;
                }
            }
        }
        ncw_close(ncid);
    }

    printf_v("  creating dst:");
    {
        int ncid_dst, ncid0;
        int ndims0;

        ncw_create(dst, NC_CLOBBER | NC_64BIT_OFFSET, &ncid_dst);

        ncw_open(src[0], NC_NOWRITE, &ncid0);
        ncw_inq_ndims(ncid0, &ndims0);

        for (i = 0; i < ndims0; ++i) {
            char dimname[NC_MAX_NAME];
            size_t dimlen;
            int dimid;
            int j;

            ncw_inq_dim(ncid0, i, dimname, &dimlen);
            for (j = 0; j < ndims; ++j) {
                if (i == dimids[0][j]) {
                    int ncid_now, dimid_now;
                    size_t dimlen_now;
                    int s;

                    for (s = 1; s < nsrc; ++s) {
                        ncw_open(src[s], NC_NOWRITE, &ncid_now);
                        ncw_inq_dimid(ncid_now, dimname, &dimid_now);
                        ncw_inq_dimlen(ncid_now, dimid_now, &dimlen_now);
                        ncw_close(ncid_now);
                        dimlen += dimlen_now;
                    }
                    break;
                }
            }
            ncw_def_dim(ncid_dst, dimname, dimlen, &dimid);
        }

        for (i = 0; i < nvars0; ++i) {
            nc_type type;
            int nvardims;
            char varname[NC_MAX_NAME];
            int dimids0[NC_MAX_DIMS], dimids_dst[NC_MAX_DIMS];
            int natts;
            int varid_dst;
            int j;

            ncw_inq_var(ncid0, i, varname, &type, &nvardims, dimids0, &natts);
            for (j = 0; j < nvardims; ++j) {
                char dimname[NC_MAX_NAME];
                size_t dimlen;

                ncw_inq_dim(ncid0, dimids0[j], dimname, &dimlen);
                ncw_inq_dimid(ncid_dst, dimname, &dimids_dst[j]);
            }

            ncw_def_var(ncid_dst, varname, type, nvardims, dimids_dst, &varid_dst);
            ncw_copy_atts(ncid0, i, ncid_dst, varid_dst);
        }

        ncw_close(ncid_dst);
        ncw_close(ncid0);
    }
    printf_v("\n");

    printf_v("  copying data:");
    size = malloc(nsrc * sizeof(size_t));
    {
        int ncid_dst;

        ncw_open(dst, NC_WRITE, &ncid_dst);

        for (i = 0; i < nvars0; ++i) {
            int j;

            for (j = 0; j < nvars; ++j)
                if (varids[j] == i)
                    break;

            if (j == nvars) {
                int ncid0;

                ncw_open(src[0], NC_NOWRITE, &ncid0);
                ncw_copy_vardata(ncid0, i, ncid_dst);
                ncw_close(ncid0);
            } else {
                int ncid;
                char varname[NC_MAX_NAME];
                int s;
                size_t size_total = 0;
                size_t nbytes;
                void* data = NULL;
                void* data_now;
                nc_type type;
                int varid_dst;

                ncw_open(src[0], NC_NOWRITE, &ncid);
                ncw_inq_varname(ncid, i, varname);
                ncw_close(ncid);

                size_total = 0;
                for (s = 0; s < nsrc; ++s) {
                    int ncid;
                    int varid;
                    int ndims;
                    int dimids[NC_MAX_DIMS];
                    size_t dimlens[NC_MAX_DIMS];

                    ncw_open(src[s], NC_NOWRITE, &ncid);
                    ncw_inq_varid(ncid, varname, &varid);
                    ncw_inq_var(ncid, varid, NULL, &type, &ndims, dimids, NULL);
                    for (j = 0; j < ndims; ++j)
                        ncw_inq_dimlen(ncid, dimids[j], &dimlens[j]);
                    ncw_close(ncid);
                    size[s] = 1;
                    for (j = 0; j < ndims; ++j)
                        size[s] *= dimlens[j];
                    size_total += size[s];
                }

                nbytes = size_total * ncw_sizeof(type);
                if (nbytes > 4294967296)
                    quit("sizeof(%s) = %zu exceeds 4GB", varname, nbytes);
                data = malloc(nbytes);
                if (data == NULL)
                    quit("malloc(): could not allocate memory for variable \"%s\" (size = %zu)", varname, size_total * ncw_sizeof(type));
                data_now = data;
                for (s = 0; s < nsrc; ++s) {
                    int ncid;
                    int varid;

                    ncw_open(src[s], NC_NOWRITE, &ncid);
                    ncw_inq_varid(ncid, varname, &varid);
                    ncw_get_var(ncid, varid, data_now);
                    ncw_close(ncid);
                    data_now = &((char*) data_now)[size[s] * ncw_sizeof(type)];
                }
                ncw_inq_varid(ncid_dst, varname, &varid_dst);
                ncw_put_var(ncid_dst, varid_dst, data);
                free(data);
            }
            printf_v(".");
        }
        ncw_close(ncid_dst);
    }
    printf_v("\n");

    /*
     * clean up
     */
    free(size);
    free(varids);
    free(dimids);
    for (i = 0; i < nsrc; ++i)
        free(src[i]);
    free(src);

    return 0;
}
Esempio n. 7
0
void das_getHE(dasystem* das)
{
    observations* obs = das->obs;
    model* m = das->m;
    ENSOBSTYPE* Hx = NULL;
    int i, e;

    das->s_mode = S_MODE_HE_f;
    if (obs->nobs == 0)
        return;

    if (das->nmem <= 0)
        das_setnmem(das);
    enkf_printf("    ensemble size = %d\n", das->nmem);
    assert(das->nmem > 0);

    distribute_iterations(0, das->nmem - 1, nprocesses, rank, "    ");

    /*
     * ensemble observation array to be filled 
     */
    assert(das->S == NULL);
    das->S = alloc2d(das->nmem, obs->nobs, sizeof(ENSOBSTYPE));
    if (das->mode == MODE_ENOI)
        Hx = calloc(obs->nobs, sizeof(ENSOBSTYPE));

    for (i = 0; i < obs->nobstypes; ++i) {
        obstype* ot = &obs->obstypes[i];
        float*** vvv = NULL;
        float** vv = NULL;
        H_fn H = NULL;
        int mvid;
        int ni, nj, nk;
        int nobs;
        int* obsids;
        char fname[MAXSTRLEN];

        enkf_printf("    %s ", ot->name);
        fflush(stdout);

        mvid = model_getvarid(m, obs->obstypes[i].varnames[0], 1);
        if (ot->issurface) {
            model_getvardims(m, mvid, &ni, &nj, NULL);
            vv = alloc2d(nj, ni, sizeof(float));
        } else {
            model_getvardims(m, mvid, &ni, &nj, &nk);
            vvv = alloc3d(nk, nj, ni, sizeof(float));
        }

        /*
         * set H
         */
        H = getH(ot->name, ot->hfunction);

        if (ot->isasync) {
            int t1 = get_tshift(ot->date_min, ot->async_tstep);
            int t2 = get_tshift(ot->date_max, ot->async_tstep);
            int t;

            for (t = t1; t <= t2; ++t) {
                enkf_printf("|");
                obs_find_bytypeandtime(obs, i, t, &nobs, &obsids);
                if (nobs == 0)
                    continue;

                /*
                 * for EnOI it is essential sometimes (e.g. in some bias
                 * correction cases) that the background is interpolated first
                 */
                if (das->mode == MODE_ENOI) {
                    if (enkf_obstype == OBSTYPE_VALUE) {
                        int success = model_getbgfname_async(m, das->bgdir, ot->varnames[0], ot->name, t, fname);

                        H(das, nobs, obsids, fname, -1, t, (ot->issurface) ? (void*) vv : (void*) vvv, Hx);
                        enkf_printf((success) ? "A" : "S");
                        fflush(stdout);
                    } else if (enkf_obstype == OBSTYPE_INNOVATION) {
                        Hx[0] = 0;
                        enkf_printf("-");
                        fflush(stdout);
                    }
                }

                if (das->mode == MODE_ENKF || !enkf_fstatsonly) {
                    for (e = my_first_iteration; e <= my_last_iteration; ++e) {
                        int success = model_getmemberfname_async(m, das->ensdir, ot->varnames[0], ot->name, e + 1, t, fname);

                        H(das, nobs, obsids, fname, e + 1, t, (ot->issurface) ? (void*) vv : (void*) vvv, das->S[e]);
                        enkf_printf((success) ? "a" : "s");
                        fflush(stdout);
                    }
                }

                free(obsids);
            }
        } else {
            obs_find_bytype(obs, i, &nobs, &obsids);
            if (nobs == 0)
                goto next;

            /*
             * for EnOI it is essential sometimes (e.g. in some bias correction
             * cases) that the background is interpolated first
             */
            if (das->mode == MODE_ENOI) {
                if (enkf_obstype == OBSTYPE_VALUE) {
                    model_getbgfname(m, das->bgdir, ot->varnames[0], fname);
                    H(das, nobs, obsids, fname, -1, INT_MAX, (ot->issurface) ? (void*) vv : (void*) vvv, Hx);
                    enkf_printf("+");
                    fflush(stdout);
                } else if (enkf_obstype == OBSTYPE_INNOVATION) {
                    Hx[0] = 0;
                    enkf_printf("-");
                    fflush(stdout);
                }
            }

            if (das->mode == MODE_ENKF || !enkf_fstatsonly) {
                for (e = my_first_iteration; e <= my_last_iteration; ++e) {
                    model_getmemberfname(m, das->ensdir, ot->varnames[0], e + 1, fname);
                    H(das, nobs, obsids, fname, e + 1, INT_MAX, (ot->issurface) ? (void*) vv : (void*) vvv, das->S[e]);
                    enkf_printf(".");
                    fflush(stdout);
                }
            }

            free(obsids);
        }

      next:

        if (ot->issurface)
            free(vv);
        else
            free(vvv);
        enkf_printf("\n");
    }                           /* for i (over obstypes) */

#if defined(MPI)
    if (das->mode == MODE_ENKF || !enkf_fstatsonly) {
#if !defined(HE_VIAFILE)
        /*
         * communicate HE via MPI
         */
        int ierror, sendcount, *recvcounts, *displs;

        recvcounts = malloc(nprocesses * sizeof(int));
        displs = malloc(nprocesses * sizeof(int));

        sendcount = my_number_of_iterations * obs->nobs;
        for (i = 0; i < nprocesses; ++i) {
            recvcounts[i] = number_of_iterations[i] * obs->nobs;
            displs[i] = first_iteration[i] * obs->nobs;
        }

        ierror = MPI_Allgatherv(das->S[my_first_iteration], sendcount, MPIENSOBSTYPE, das->S[0], recvcounts, displs, MPIENSOBSTYPE, MPI_COMM_WORLD);
        assert(ierror == MPI_SUCCESS);

        free(recvcounts);
        free(displs);
#else
        /*
         * communicate HE via file
         */
        {
            int ncid;
            int varid;
            size_t start[2], count[2];

            if (rank == 0) {
                int dimids[2];

                ncw_create(FNAME_HE, NC_CLOBBER | NETCDF_FORMAT, &ncid);
                ncw_def_dim(FNAME_HE, ncid, "m", das->nmem, &dimids[0]);
                ncw_def_dim(FNAME_HE, ncid, "p", obs->nobs, &dimids[1]);
                ncw_def_var(FNAME_HE, ncid, "HE", NC_FLOAT, 2, dimids, &varid);
                ncw_close(FNAME_HE, ncid);
            }
            MPI_Barrier(MPI_COMM_WORLD);

            ncw_open(FNAME_HE, NC_WRITE, &ncid);
            ncw_inq_varid(FNAME_HE, ncid, "HE", &varid);
            start[0] = my_first_iteration;
            start[1] = 0;
            count[0] = my_last_iteration - my_first_iteration + 1;
            count[1] = obs->nobs;
            ncw_put_vara_float(FNAME_HE, ncid, varid, start, count, das->S[my_first_iteration]);
            ncw_close(FNAME_HE, ncid);
            MPI_Barrier(MPI_COMM_WORLD);

            ncw_open(FNAME_HE, NC_NOWRITE, &ncid);
            ncw_inq_varid(FNAME_HE, ncid, "HE", &varid);
            ncw_get_var_float(FNAME_HE, ncid, varid, das->S[0]);
            ncw_close(FNAME_HE, ncid);
        }
#endif
    }
#endif

    if (das->mode == MODE_ENOI) {
        /*
         * subtract ensemble mean; add background
         */
        if (!enkf_fstatsonly) {
            double* ensmean = calloc(obs->nobs, sizeof(double));

            for (e = 0; e < das->nmem; ++e) {
                ENSOBSTYPE* Se = das->S[e];

                for (i = 0; i < obs->nobs; ++i)
                    ensmean[i] += Se[i];
            }
            for (i = 0; i < obs->nobs; ++i)
                ensmean[i] /= (double) das->nmem;

            for (e = 0; e < das->nmem; ++e) {
                ENSOBSTYPE* Se = das->S[e];

                for (i = 0; i < obs->nobs; ++i)
                    Se[i] += Hx[i] - ensmean[i];
            }

            free(ensmean);
        } else {
            for (e = 0; e < das->nmem; ++e) {
                ENSOBSTYPE* Se = das->S[e];

                for (i = 0; i < obs->nobs; ++i)
                    Se[i] = Hx[i];
            }
        }
    }

    if (das->mode == MODE_ENOI)
        free(Hx);
}