File: triangle2d.c

package info (click to toggle)
adios 1.13.1-31
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 23,692 kB
  • sloc: ansic: 133,236; f90: 8,791; sh: 7,779; python: 7,648; xml: 3,793; makefile: 2,996; cpp: 2,340; java: 626; sed: 16; perl: 8
file content (112 lines) | stat: -rw-r--r-- 3,496 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
/* 
 * ADIOS is freely available under the terms of the BSD license described
 * in the COPYING file in the top level directory of this source distribution.
 *
 * Copyright (c) 2008 - 2009.  UT-BATTELLE, LLC. All rights reserved.
 */

/* ADIOS C Example: write variables along with an unstructured mesh. 
 */
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "mpi.h"
#include "adios.h"




int main (int argc, char ** argv ) 
{
    int         rank, size, i, y;
    MPI_Comm    comm = MPI_COMM_WORLD;

    int         npoints, num_cells;
    float       *points;
    int         *cells;
    double      *N; // node centered variable
    double      *C; // cell centered variable
	

    /* ADIOS variables declarations for matching gwrite_temperature.ch */
    uint64_t    adios_groupsize, adios_totalsize;
    int64_t     adios_handle;

    MPI_Init (&argc, &argv);
    MPI_Comm_rank (comm, &rank);
    MPI_Comm_size (comm, &size);

    int n = size; // number of points in y direction
    npoints = n*(size+1);       // global number of points
    num_cells = 2*(n-1)*size;   // global number of cells

    // local mesh + data
    points = malloc (sizeof(float) * 2*n*2);  // 2n points with 2 (x-y) coords
    cells = malloc (sizeof(int) * 2*(n-1)*3); // 2(n-1) triangles
    N = malloc (sizeof(double) * 2*n);           // 2n data points on the points
    C = malloc (sizeof(double) * 2*(n-1));       // 2(n-1) data points on the cells

    // generate local data 
    // points (rank,0), (rank,1),...,(rank,n-1), 
    //        (rank+1,0) ... (rank+1,n-1)
    int lp = 2*n; // number of points in this local array
    int op = n*rank; // offset in global points array (next process overwrites our second half)
    for (y=0; y<n; y++) {
        points[2*y] = rank;
        points[2*y+1] = y+rank*0.1;
    }
    // second half are ghost cells with next process
    for (y=0; y<n; y++) {
        points[2*n+2*y] = rank+1;
        points[2*n+2*y+1] = y+(rank+1)*0.1;
    }

    // cells (0--n--n+1) (0--n+1--1) (1--n+1--n+2) (1--n+2--2)... (n-2, 2n-1, n-1) in local terms
    // 0 point local is rank*n in global = op
    int lc = 2*(n-1);
    int oc = lc*rank;
    for (i=0; i<n-1; i++) {
        // two triangles define one rectangle 
        int p=6*i;
        cells[p+0] = op+i; cells[p+1] = op+n+i;   cells[p+2] = op+n+i+1;
        cells[p+3] = op+i; cells[p+4] = op+n+i+1; cells[p+5] = op+i+1;
    }

    for (i=0; i<lp; i++)
        N[i] = rank*n+i;

    for (i=0; i<lc; i++)
        C[i] = rank*lc+i;

    adios_init ("triangle2d.xml", comm);

    adios_open (&adios_handle, "tri", "triangle2d.bp", "w", comm);

    adios_groupsize = 7*sizeof(int) \
	+ sizeof(double) * (lp) \
	+ sizeof(double) * (lc) \
	+ sizeof(float) * (lp) * (2) \
	+ sizeof(int) * (lc) * (3);

    adios_group_size (adios_handle, adios_groupsize, &adios_totalsize);
    adios_write (adios_handle, "nproc", &size);
    adios_write (adios_handle, "npoints", &npoints);
    adios_write (adios_handle, "num_cells", &num_cells);
    adios_write (adios_handle, "lp", &lp);
    adios_write (adios_handle, "op", &op);
    adios_write (adios_handle, "lc", &lc);
    adios_write (adios_handle, "oc", &oc);
    adios_write (adios_handle, "cells", cells);
    adios_write (adios_handle, "points", points);
    adios_write (adios_handle, "N", N);
    adios_write (adios_handle, "C", C);

    adios_close (adios_handle);

    MPI_Barrier (comm);

    adios_finalize (rank);

    MPI_Finalize ();
    return 0;
}