~njansson/dolfin/hpc

« back to all changes in this revision

Viewing changes to src/bench/fem/passembly/main.cpp

  • Committer: Johannes Ring
  • Date: 2008-03-05 22:43:06 UTC
  • Revision ID: johannr@simula.no-20080305224306-2npsdyhfdpl2esji
The BIG commit!

Show diffs side-by-side

added added

removed removed

Lines of Context:
1
 
// Copyright (C) 2007 Garth N. Wells.
2
 
// Licensed under the GNU LGPL Version 2.1.
3
 
//
4
 
// Modified by Magnus Vikstrøm 
5
 
// First added:  2007-10-30
6
 
// Last changed: 2008-01-09
7
 
//
8
 
// This file is used for testing distribution of the mesh using MPI
9
 
 
10
 
#include <dolfin.h>
11
 
#include "Poisson.h"
12
 
#include "Nonlinear2D.h"
13
 
#include <iostream>
14
 
#include <fstream>
15
 
#include "getopts.h"
16
 
#include "assemble.h"
17
 
 
18
 
using namespace dolfin;
19
 
 
20
 
real timer(Mesh& mesh, int num_iterations, Form& a)
21
 
{
22
 
  dolfin::cout << "Assembling with sequential assembler." << dolfin::endl;
23
 
  Matrix A;
24
 
  Assembler assembler(mesh);
25
 
 
26
 
  assembler.assemble(A, a, true);
27
 
  tic();
28
 
  for(int i=0; i<num_iterations; ++i)
29
 
    assembler.assemble(A, a, false);
30
 
  return toc()/static_cast<real>(num_iterations);
31
 
}
32
 
 
33
 
real p_timer(Mesh& mesh, MeshFunction<dolfin::uint>& partitions, int num_iterations, Form& a)
34
 
{
35
 
  dolfin::cout << "Assembling with parallel assembler." << dolfin::endl;
36
 
  Matrix B;
37
 
  pAssembler passembler(mesh, partitions);
38
 
 
39
 
  passembler.assemble(B, a, true);
40
 
  tic();
41
 
  for(int i=0; i<num_iterations; ++i)
42
 
    passembler.assemble(B, a, false);
43
 
  return toc()/static_cast<real>(num_iterations);
44
 
}
45
 
 
46
 
int main(int argc, char* argv[])
47
 
{
48
 
  std::string meshfile, partitionsfile, resultfile;
49
 
  std::string assembler = "parallel";
50
 
  int cells = 400;
51
 
  int cells_3D = 40;
52
 
  int num_iterations = 1;
53
 
  int num_part = dolfin::MPI::numProcesses();
54
 
  int debug = -1;
55
 
  bool check = false;
56
 
  bool sequential = false;
57
 
  std::string testtype = "Poisson3D";
58
 
  bool petsc_info = false;
59
 
 
60
 
  Options listOpts;
61
 
  int switchInt;
62
 
  listOpts.addOption("", "mesh", "Mesh File", true);
63
 
  listOpts.addOption("", "partitions", "Mesh partitions file", true);
64
 
  listOpts.addOption("", "sequential", "Run with sequential assembler", false);
65
 
  listOpts.addOption("", "cells", "Number of cells", true);
66
 
  listOpts.addOption("", "num_part", "Number of partitions", true);
67
 
  listOpts.addOption("", "resultfile", "File to save results in", true);
68
 
  listOpts.addOption("", "num_iterations", "Number of times to assemble", true);
69
 
  listOpts.addOption("", "check", "Verify assembly result", false);
70
 
  listOpts.addOption("", "debug", "Prints debugging info", false);
71
 
  listOpts.addOption("", "testtype", "Type of test: Poisson2D, Poisson3D, Nonlinear2D", true);
72
 
  listOpts.addOption("", "petsc_info", "Print info from PETSc", false);
73
 
 
74
 
  if (listOpts.parse(argc, argv))
75
 
    while ((switchInt = listOpts.cycle()) >= 0)
76
 
    {
77
 
      switch(switchInt)
78
 
      {
79
 
        case 0:
80
 
          meshfile = listOpts.getArgs(switchInt);
81
 
          break;
82
 
        case 1:
83
 
          partitionsfile = listOpts.getArgs(switchInt);
84
 
          break;
85
 
        case 2:
86
 
          sequential = true;
87
 
          break;
88
 
        case 3:
89
 
          cells = atoi(listOpts.getArgs(switchInt).c_str());
90
 
          cells_3D = atoi(listOpts.getArgs(switchInt).c_str());
91
 
          break;
92
 
        case 4:
93
 
          num_part = atoi(listOpts.getArgs(switchInt).c_str());
94
 
          break;
95
 
        case 5:
96
 
          resultfile = listOpts.getArgs(switchInt);
97
 
          break;
98
 
        case 6:
99
 
          num_iterations = atoi(listOpts.getArgs(switchInt).c_str());
100
 
          break;
101
 
        case 7:
102
 
          check = true;
103
 
          break;
104
 
        case 8:
105
 
          debug = 1;
106
 
          break;
107
 
        case 9:
108
 
          testtype = listOpts.getArgs(switchInt);
109
 
          break;
110
 
        case 10:
111
 
          petsc_info = true;
112
 
          break;
113
 
        default:
114
 
          break;
115
 
      }
116
 
    }
117
 
 
118
 
  Mesh mesh;
119
 
  MeshFunction<dolfin::uint>* partitions;
120
 
  dolfin_set("debug level", debug);
121
 
  if(petsc_info)
122
 
  {
123
 
    char** init_argv = new char*[2];
124
 
    init_argv[0] = argv[0];
125
 
    sprintf(init_argv[1], "%s", "-info");
126
 
    dolfin::cout << "dolfin_init" << dolfin::endl;
127
 
    dolfin_init(2, init_argv);
128
 
  }
129
 
  if(meshfile != "")
130
 
  {
131
 
    dolfin::cout << "Reading mesh from file: " << meshfile << dolfin::endl;
132
 
    mesh = Mesh(meshfile);
133
 
  }
134
 
  else
135
 
  {
136
 
    if(testtype == "Nonlinear2D" || testtype == "Poisson2D")
137
 
    {
138
 
      printf("Creating UnitSquare(%d, %d)\n",  cells, cells);
139
 
      mesh = UnitSquare(cells, cells);
140
 
    }
141
 
    else
142
 
    {
143
 
      printf("Creating UnitCube(%d, %d, %d)\n", cells_3D, cells_3D, cells_3D);
144
 
      mesh = UnitCube(cells_3D, cells_3D, cells_3D);
145
 
    }
146
 
  }
147
 
  real time = 0;
148
 
  if(sequential)
149
 
  {
150
 
    Form* a;
151
 
    if(testtype == "Poisson3D")
152
 
      a = new PoissonBilinearForm();
153
 
    else if(testtype == "Nonlinear2D")
154
 
    {
155
 
      Function w(mesh, 1.0);
156
 
      a = new Nonlinear2DBilinearForm(w);
157
 
    }
158
 
    else
159
 
      a = new PoissonBilinearForm();
160
 
 
161
 
    dolfin::cout << "Running test " << testtype << dolfin::endl;
162
 
    dolfin::cout << "Assembling with sequential assembler." << dolfin::endl;
163
 
    printf("Number of iteration(s): %d\n", num_iterations);
164
 
    time = timer(mesh, num_iterations, *a);
165
 
  }
166
 
  else
167
 
  {
168
 
    if(partitionsfile != "")
169
 
    {
170
 
      dolfin::cout << "Reading partitions from file: " << partitionsfile << dolfin::endl;
171
 
      partitions = new MeshFunction<dolfin::uint>(mesh, partitionsfile);
172
 
    }
173
 
    else
174
 
    {
175
 
      dolfin::cout << "Partitioning mesh into: " << num_part << " partitions" << dolfin::endl;
176
 
      partitions = new MeshFunction<dolfin::uint>(mesh);
177
 
      mesh.partition(*partitions, num_part);
178
 
    }
179
 
 
180
 
    Form* a;
181
 
    if(testtype == "Poisson3D")
182
 
      a = new PoissonBilinearForm();
183
 
    else if(testtype == "Nonlinear2D")
184
 
    {
185
 
      Function w(mesh, 1.0);
186
 
      a = new Nonlinear2DBilinearForm(w);
187
 
    }
188
 
    else
189
 
      a = new PoissonBilinearForm();
190
 
 
191
 
    dolfin::cout << "Running test " << testtype << dolfin::endl;
192
 
    dolfin::cout << "Assembling with parallel assembler." << dolfin::endl;
193
 
    printf("Number of iteration(s): %d Number of partitions: %d\n", num_iterations, num_part);
194
 
    time = p_timer(mesh, *partitions, num_iterations, *a);
195
 
  }
196
 
  if(resultfile != "")
197
 
  {
198
 
    dolfin::cout << "Appending results to " << resultfile << dolfin::endl;
199
 
    std::ofstream outfile(resultfile.c_str(), std::ofstream::app);
200
 
    outfile << dolfin::MPI::numProcesses() << " " << time << std::endl;
201
 
    outfile.close();
202
 
  }
203
 
  else
204
 
  {
205
 
    printf("Average assemble time: %.3e\n", time);
206
 
  }
207
 
  if(check)
208
 
  {
209
 
    dolfin::cout << "Not implemented" << dolfin::endl;
210
 
    //check_assembly(mesh, *partitions);
211
 
  }
212
 
  return 0;
213
 
}