File: test_nodes2.c

package info (click to toggle)
p4est 2.3.6-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 8,536 kB
  • sloc: ansic: 87,528; makefile: 855; sh: 635; perl: 272; python: 226; awk: 40; javascript: 23
file content (103 lines) | stat: -rw-r--r-- 2,974 bytes parent folder | download | duplicates (2)
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
/*
  This file is part of p4est.
  p4est is a C library to manage a collection (a forest) of multiple
  connected adaptive quadtrees or octrees in parallel.

  Copyright (C) 2010 The University of Texas System
  Additional copyright (C) 2011 individual authors
  Written by Carsten Burstedde, Lucas C. Wilcox, and Tobin Isaac

  p4est is free software; you can redistribute it and/or modify
  it under the terms of the GNU General Public License as published by
  the Free Software Foundation; either version 2 of the License, or
  (at your option) any later version.

  p4est is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  GNU General Public License for more details.

  You should have received a copy of the GNU General Public License
  along with p4est; if not, write to the Free Software Foundation, Inc.,
  51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/

#ifndef P4_TO_P8
#include <p4est.h>
#include <p4est_extended.h>
#include <p4est_nodes.h>
#else
#include <p8est.h>
#include <p8est_extended.h>
#include <p8est_nodes.h>
#endif

static int
refine_fn (p4est_t * p4est, p4est_topidx_t which_tree,
           p4est_quadrant_t * quadrant)
{
  if (quadrant->level >= P4EST_QMAXLEVEL) {
    return 0;
  }
  if (quadrant->x == 0 && quadrant->y == 0 &&
#ifdef P4_TO_P8
      quadrant->z == 0 &&
#endif
      1) {
    return 1;
  }

  return 0;
}

int
main (int argc, char **argv)
{
  int                 mpiret;
  sc_MPI_Comm         mpicomm;
  p4est_t            *p4est;
  p4est_connectivity_t *connectivity;
  p4est_ghost_t      *ghost;
  p4est_nodes_t      *nodes1, *nodes2;

  mpiret = sc_MPI_Init (&argc, &argv);
  SC_CHECK_MPI (mpiret);
  mpicomm = sc_MPI_COMM_WORLD;

  sc_init (mpicomm, 1, 1, NULL, SC_LP_DEFAULT);

  /* create connectivity and forest structures */
#ifndef P4_TO_P8
  connectivity = p4est_connectivity_new_unitsquare ();
#else
  connectivity = p8est_connectivity_new_unitcube ();
#endif
  /* specifying min_quadrants is not partition independent but fun */
  p4est = p4est_new_ext (mpicomm, connectivity, 15, 0, 0, 0, NULL, NULL);

  /* refine and balance to make the number of elements interesting */
  p4est_refine (p4est, 1, refine_fn, NULL);
  p4est_balance (p4est, P4EST_CONNECT_FULL, NULL);

  /* create the ghost structure */
  ghost = p4est_ghost_new (p4est, P4EST_CONNECT_FULL);

  /* create the nodes structure */
  P4EST_GLOBAL_INFO ("Making nodes with ghosts\n");
  nodes1 = p4est_nodes_new (p4est, ghost);
  P4EST_GLOBAL_INFO ("Making nodes without ghosts\n");
  nodes2 = p4est_nodes_new (p4est, NULL);

  /* clean up and exit */
  p4est_nodes_destroy (nodes1);
  p4est_nodes_destroy (nodes2);
  p4est_ghost_destroy (ghost);
  p4est_destroy (p4est);
  p4est_connectivity_destroy (connectivity);
  sc_finalize ();

  mpiret = sc_MPI_Finalize ();
  SC_CHECK_MPI (mpiret);

  return 0;
}