Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
advectest.C
Go to the documentation of this file.
1 /* *******************************************************************
2  * Rocstar Simulation Suite *
3  * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4  * *
5  * Illinois Rocstar LLC *
6  * Champaign, IL *
7  * www.illinoisrocstar.com *
8  * sales@illinoisrocstar.com *
9  * *
10  * License: See LICENSE file in top level of distribution package or *
11  * http://opensource.org/licenses/NCSA *
12  *********************************************************************/
13 /* *******************************************************************
14  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15  * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16  * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17  * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18  * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20  * ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21  * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22  *********************************************************************/
23 // $Id: advectest.C,v 1.6 2008/12/06 08:45:28 mtcampbe Exp $
24 
25 #include "roccom.h"
26 
27 #include <cstdio>
28 #include <iostream>
29 #include <fstream>
30 #include <algorithm>
31 #include <cstring>
32 #include <string>
33 #include <cstdlib>
34 #include <cmath>
35 #include <cassert>
36 #include <sstream>
37 #include "roccom_assertion.h"
38 #include "../Rocblas/include/Rocblas.h"
39 #include "../Rocsurf/include/Rocsurf.h"
40 
41 #include "../Rocsurf/test/IM_Reader.h"
42 #include "PointPropagate.h"
43 
44 #ifdef MESH_ADAPT
45 #include "AdaptCOMWindow.h"
46 #endif
47 
48 
49 using namespace std;
50 using namespace PROP;
51 
53 COM_EXTERN_MODULE( Rocmap);
57 
58 void load_modules() {
60  COM_LOAD_MODULE_STATIC_DYNAMIC(Rocmap, "MAP");
64 }
65 
66 static int rank = 0;
67 
68 void print_usage( int argc, char *argv[]) {
69  if ( argc <= 2) {
70  cout << "Usage: " << argv[0] << " <surffile> <controlfile>" << std::endl;
71 
72  exit(-1);
73  }
74 }
75 
77  Control_parameter() : speed(0), timestep(0), steps(0), start(0), interval(0),
78  remesh_interval(10), test("rotate"),
79  adapt_iter(2), do_redist(0), do_collapse(1), do_flip(1),
80  do_split(1), collapse_ratio(0), split_angle(0), refine(1) {}
81 
82  string method;
83  string wavefrontal;
84  string normaldif;
85  string eigthres;
86  string courant;
87  string fangle;
88  string sangle;
89  string smoother;
90  string rediter;
91  string weight;
92  string conserv;
93 
94  double speed;
95  double timestep;
96  int steps;
97  int start;
98  int interval;
100 
101  string test;
102  string verbose;
103 
104  // The following options are for mesh adaptation
108  int do_flip;
109  int do_split;
111  double split_angle;
112  double refine;
113 };
114 
115 void read_control_file( const char *fname, Control_parameter &cp) {
116  /* Sample control file:
117  * method: fo # method: fo and mp
118  * wavefrontal: 1 # wavefrontal expansion
119  * normaldif: 1 # normal diffusion
120  * eigthres: 1.e-4 # threshold for null space: 0..1 (1.e-4 default)
121  * courant: 0.5 # courant constant
122  * fangle: 25 # feature edge angle: between 0 and 180
123  * sangle: 65 # feature edge angle: between 0 and 180
124  * smoother: angle # type of mesh-smoothing algorithm
125  * rediter: 1 # Number of iterations for vertex redistribution
126  * conserv: 0 # Whether or not to conserve volume
127  * speed: 0.1 # Speed
128  * timestep: 0.001 # time step
129  * steps: 100 # number of time steps
130  * interval: 10 # output intervals
131  * remesh_interval: 10 # remesh intervals
132  * weight: 3 # angle weighted
133  * verbose: 1 # verbose level
134  * test: rotate # test case
135  *
136  * adapt_iter: 2
137  * do_flip: 1
138  * do_split: 1
139  * do_redist: 0
140  * do_collapse: 1
141  * collapse_ratio: 0.2
142  * split_angle: 2.618 # (150 degrees)
143  */
144  ifstream is(fname); COM_assertion_msg( is, "File does not exist");
145 
146  while ( !is.eof()) {
147  char buf[255];
148  is.get( buf, 255, ':');
149  if ( buf[0] == '\0') { is.getline( buf, 255); continue; }
150 
151  istringstream istr(buf);
152  string keyword; istr >> keyword;
153  is.getline( buf, 255, ':');
154 
155  if ( keyword == "method")
156  is >> cp.method;
157  else if ( keyword == "wavefrontal")
158  is >> cp.wavefrontal;
159  else if ( keyword == "normaldif")
160  is >> cp.normaldif;
161  else if ( keyword == "eigthres")
162  is >> cp.eigthres;
163  else if ( keyword == "courant")
164  is >> cp.courant;
165  else if ( keyword == "fangle")
166  is >> cp.fangle;
167  else if ( keyword == "sangle")
168  is >> cp.sangle;
169  else if ( keyword == "smoother")
170  is >> cp.smoother;
171  else if ( keyword == "rediter")
172  is >> cp.rediter;
173  else if ( keyword == "conserv")
174  is >> cp.conserv;
175  else if ( keyword == "speed")
176  is >> cp.speed;
177  else if ( keyword == "timestep")
178  is >> cp.timestep;
179  else if ( keyword == "steps")
180  is >> cp.steps;
181  else if ( keyword == "start")
182  is >> cp.start;
183  else if ( keyword == "interval")
184  is >> cp.interval;
185  else if ( keyword == "remesh_interval")
186  is >> cp.remesh_interval;
187  else if ( keyword == "weight")
188  is >> cp.weight;
189  else if ( keyword == "test")
190  is >> cp.test;
191  else if ( keyword == "verbose")
192  is >> cp.verbose;
193  else if ( keyword == "adapt_iter")
194  is >> cp.adapt_iter;
195  else if ( keyword == "do_redist")
196  is >> cp.do_redist;
197  else if ( keyword == "do_contract" || keyword == "do_collapse")
198  is >> cp.do_collapse;
199  else if ( keyword == "do_flip")
200  is >> cp.do_flip;
201  else if ( keyword == "do_split")
202  is >> cp.do_split;
203  else if ( keyword == "refine")
204  is >> cp.refine;
205  else if ( keyword == "collapse_ratio")
206  is >> cp.collapse_ratio;
207  else if ( keyword == "split_angle")
208  is >> cp.split_angle;
209  else
210  std::cerr << "Unknow keyword " << keyword << std::endl;
211  is.getline( buf, 255);
212  }
213 
214  if ( rank==0) std::cout << " speed is " << cp.speed << std::endl;
215 }
216 
217 void init_parameters( const Control_parameter &cntr_param) {
218  int PROP_set_option = COM_get_function_handle( "PROP.set_option");
219 
220  if ( !cntr_param.method.empty()) {
221  COM_call_function( PROP_set_option, "method", cntr_param.method.c_str());
222 
223  if ( rank==0) std::cout << "Set propagation method to " << cntr_param.method << std::endl;
224  }
225 
226  if ( !cntr_param.wavefrontal.empty()) {
227  COM_call_function( PROP_set_option, "wavefrontal", cntr_param.wavefrontal.c_str());
228  if ( rank==0) std::cout << "Set wavefrontal to " << cntr_param.wavefrontal << std::endl;
229  }
230 
231  if ( !cntr_param.normaldif.empty()) {
232  COM_call_function( PROP_set_option, "normaldif", cntr_param.normaldif.c_str());
233  if ( rank==0) std::cout << "Set normaldif to " << cntr_param.normaldif << std::endl;
234  }
235 
236  if ( !cntr_param.eigthres.empty()) {
237  COM_call_function( PROP_set_option, "eigthres", cntr_param.eigthres.c_str());
238  if ( rank==0) std::cout << "Set eigthres to " << cntr_param.eigthres << std::endl;
239  }
240 
241  if ( !cntr_param.courant.empty()) {
242  COM_call_function( PROP_set_option, "courant", cntr_param.courant.c_str());
243  if ( rank==0) std::cout << "Set courant constant to " << cntr_param.courant << std::endl;
244  }
245 
246  if ( !cntr_param.fangle.empty()) {
247  COM_call_function( PROP_set_option, "fangle", cntr_param.fangle.c_str());
248  if ( rank==0) std::cout << "Set weak-feature angle to " << cntr_param.fangle << std::endl;
249  }
250  if ( !cntr_param.sangle.empty()) {
251  COM_call_function( PROP_set_option, "sangle", cntr_param.sangle.c_str());
252  if ( rank==0) std::cout << "Set strong-feature angle to " << cntr_param.sangle << std::endl;
253  }
254 
255  if ( !cntr_param.smoother.empty()) {
256  COM_call_function( PROP_set_option, "smoother", cntr_param.smoother.c_str());
257  if ( rank==0) std::cout << "Set smoother to " << cntr_param.smoother << std::endl;
258  }
259 
260  if ( !cntr_param.rediter.empty()) {
261  COM_call_function( PROP_set_option, "rediter", cntr_param.rediter.c_str());
262  if ( rank==0) std::cout << "Set rediter to " << cntr_param.rediter << std::endl;
263  }
264 
265  if ( !cntr_param.conserv.empty()) {
266  COM_call_function( PROP_set_option, "conserv", cntr_param.conserv.c_str());
267  if ( rank==0) std::cout << "Set conserv to " << cntr_param.conserv << std::endl;
268  }
269 
270  if ( !cntr_param.verbose.empty()) {
271  COM_call_function( PROP_set_option, "verbose", cntr_param.verbose.c_str());
272  if ( rank==0) std::cout << "Set verbose level to " << cntr_param.verbose << std::endl;
273  }
274 
275  if ( !cntr_param.weight.empty()) {
276  COM_call_function( PROP_set_option, "weight", cntr_param.weight.c_str());
277  if ( rank==0) std::cout << "Set weight to " << cntr_param.weight << std::endl;
278  }
279 
280 }
281 
282 // Read in a surface pmesh, and return its window name.
283 std::string read_in_mesh ( const char *fname) {
284  if ( rank==0) cout << "Reading surface mesh file \"" << fname << '"' << endl;
285 
286  std::string fname_str(fname);
287 
288  std::string::size_type pos = fname_str.find_first_of( ".");
289  const string wname = fname_str.substr( 0, pos);
290 
291  if ( rank==0) cout << "Creating window \"" << wname << '"' << endl;
292 
293  IM_Reader im_reader;
294  int npanes = im_reader.read_winmesh( fname, wname, false);
295  COM_assertion_msg( npanes>=0, "Failed to read in mesh file. File empty?");
296 
297  return wname;
298 }
299 
300 void rescale_object( std::string &wname, double alpha,
301  const SURF::Vector_3<double> &origin) {
302  // Rescale and translate the object.
303  COM::Window *win=COM_get_roccom()->get_window_object(wname);
304  COM::Attribute *nc=win->attribute("nc");
305 
306  // Rescale the radius to 0.15
307  Rocblas::mul_scalar( nc, &alpha, nc);
308 
309  // Translate the origin to (0.5,0.75, 0.5);
310  COM::Attribute *x=win->attribute("1-nc");
311  COM::Attribute *y=win->attribute("2-nc");
312  COM::Attribute *z=win->attribute("3-nc");
313 
314  Rocblas::add_scalar( x, &origin[0], x);
315  Rocblas::add_scalar( y, &origin[1], y);
316  Rocblas::add_scalar( z, &origin[2], z);
317 }
318 
319 template <class T>
320 T square(T t) { return t*t; }
321 
322 
323 void init_attributes( const string &wname,
324  const Control_parameter &cntr_param) {
325 
326  COM_new_attribute((wname+".disps_total").c_str(), 'n', COM_DOUBLE, 3, "");
327  COM_new_attribute((wname+".cnstr_types_nodal").c_str(), 'n', COM_INT, 1, "");
328  COM_new_attribute((wname+".cnstr_types_facial").c_str(), 'e', COM_INT, 1, "");
329  COM_new_attribute((wname+".cnstr_types_panel").c_str(), 'p', COM_INT, 1, "");
330  COM_set_size((wname+".cnstr_types_panel").c_str(), 0, 1);
331 
332  COM_new_attribute((wname+".qvels").c_str(), 'e', COM_DOUBLE, 9, "m/s");
333  COM_new_attribute((wname+".vvels").c_str(), 'n', COM_DOUBLE, 3, "m/s");
334  COM_new_attribute((wname+".disps").c_str(), 'n', COM_DOUBLE, 3, "m");
335  COM_new_attribute((wname+".faceareas").c_str(), 'e', COM_DOUBLE, 1, "");
336 
337  // Ridges of each pane.
338  COM_new_attribute((wname+".ridges").c_str(), 'p', COM_INT, 2, "");
339  COM_set_size((wname+".ridges").c_str(), 0, 0);
340 
341  // COM_new_attribute((wname+".disps_novis").c_str(), 'n', COM_DOUBLE, 3, "");
342  // COM_new_attribute((wname+".facenormals").c_str(), 'e', COM_DOUBLE, 3, "");
343  // COM_new_attribute((wname+".facecenters").c_str(), 'e', COM_DOUBLE, 3, "");
344  // COM_new_attribute((wname+".faceheights").c_str(), 'e', COM_DOUBLE, 1, "");
345 
346  // Attribute for storing the number of eigenvalues for each node.
347  // COM_new_attribute((wname+".lambdas").c_str(), 'n', COM_DOUBLE, 3, "");
348  // COM_new_attribute((wname+".eigvecs").c_str(), 'n', COM_DOUBLE, 9, "");
349  COM_new_attribute((wname+".tangranks").c_str(), 'n', COM_CHAR, 1, "");
350  COM_new_attribute((wname+".scales").c_str(), 'n', COM_DOUBLE, 1, "");
351 
352  COM_resize_array( (wname+".atts").c_str());
353  COM_window_init_done( wname.c_str());
354 }
355 
356 // return the volume of the domain
357 double output_solution( const string &wname, const char *timelevel,
358  double ref=0.) {
359  static int OUT_write = 0, hdl;
360 
361  if ( OUT_write==0) {
362  OUT_write = COM_get_function_handle( "OUT.write_attribute");
363  hdl = COM_get_attribute_handle( (wname+".all").c_str());
364 
365  int OUT_set = COM_get_function_handle( "OUT.set_option");
366  COM_call_function( OUT_set, "format", "HDF");
367  }
368 
369  if ( timelevel) {
370  std::string fname = wname+"_"+timelevel;
371 
372  if ( !COMMPI_Initialized()) fname.append(".hdf");
373  else fname.append("_");
374 
375  COM_call_function( OUT_write, (char*)fname.c_str(),
376  &hdl, (char*)wname.c_str(), timelevel);
377  }
378 
379  static COM::Attribute *mesh=NULL;
380  if ( mesh==NULL)
381  mesh = COM_get_roccom()->get_window_object( wname)->attribute( "mesh");
382 
383  double vol;
384  SURF::Rocsurf::compute_signed_volumes( mesh, &vol);
385 
386  std::cout << "Volume of object is " << vol << std::endl;
387  if ( ref!=0) {
388  std::cout << "Relative error in volume is "
389  << std::abs((vol-ref)/ref) << std::endl;
390  }
391 
392  return vol;
393 }
394 
395 double compute_area( const string &wname) {
396  static int SURF_area = 0, BLAS_mul, BLAS_sum_scalar, area_hdl;
397 
398  if ( SURF_area==0) {
399  SURF_area = COM_get_function_handle( "SURF.compute_element_areas");
400  BLAS_mul = COM_get_function_handle( "BLAS.mul");
401  BLAS_sum_scalar = COM_get_function_handle( "BLAS.sum_scalar_MPI");
402 
403  area_hdl = COM_get_attribute_handle( (wname+".faceareas").c_str());
404  }
405 
406  COM_call_function( SURF_area, &area_hdl);
407 
408  double sum;
409  MPI_Comm comm=MPI_COMM_WORLD;
410  COM_call_function( BLAS_sum_scalar, &area_hdl, &sum, &comm);
411  return sum;
412 }
413 
414 double compute_volume( const string &wname) {
415  static int SURF_vol = 0, mesh_hdl;
416 
417  if ( SURF_vol==0) {
418  SURF_vol = COM_get_function_handle( "SURF.compute_signed_volumes");
419 
420  mesh_hdl = COM_get_attribute_handle( (wname+".mesh").c_str());
421  }
422 
423  double vol;
424  COM_call_function( SURF_vol, &mesh_hdl, &vol);
425 
426  return vol;
427 }
428 
429 int main(int argc, char *argv[]) {
430  COM_init( &argc, &argv);
431  load_modules();
432  print_usage( argc, argv);
433 
435 
436  // Read in mesh file.
437  string wname = read_in_mesh( argv[1]);
438 
439  // Read in control parameters
440  Control_parameter cntr_param;
441  read_control_file( argc>2?argv[2]:"control.txt", cntr_param);
442 
443  // Initialize the attributes and output the initial solution.
444  init_attributes( wname, cntr_param);
445 
446  // Attributes
447  int pmesh = COM_get_attribute_handle( (wname+".pmesh").c_str());
448  int nc = COM_get_attribute_handle( (wname+".nc").c_str());
449  int vvels = COM_get_attribute_handle( (wname+".vvels").c_str());
450  int qvels = COM_get_attribute_handle( (wname+".qvels").c_str());
451  int disps = COM_get_attribute_handle( (wname+".disps").c_str());
452  int disps_total = COM_get_attribute_handle( (wname+".disps_total").c_str());
453 
454  // Funcitions
455  int PROP_init = COM_get_function_handle( "PROP.initialize");
456  int PROP_propagate = COM_get_function_handle( "PROP.propagate");
457  int BLAS_add = COM_get_function_handle( "BLAS.add");
458 
459  if ( rank==0) cout << "Propagating interface..." << endl;
460  // Initialize control parameters
461  COM_call_function( PROP_init, &pmesh);
462  init_parameters( cntr_param);
463 
464  Translate_speed trans_spd( cntr_param.speed);
465  Rotate_speed rotate_spd;
466  Vortex_flow vortex_spd(cntr_param.speed);
467  LeVeque_flow leveque_spd(cntr_param.speed);
468 
469  bool use_quadpoints = false;
470  Speed *spd=NULL;
471  std::cerr << "Initializing for test case " << cntr_param.test << std::endl;
472  if ( cntr_param.test == "trans") {
473  spd = &trans_spd;
474  }
475  else if ( cntr_param.test == "rotate") {
476  // use_quadpoints = true;
477  spd = &rotate_spd;
478  rescale_object( wname, 1, SURF::Vector_3<double>(5./3., 0, 0));
479  }
480  else if ( cntr_param.test == "vortexcube") {
481  spd = &vortex_spd;
482  if ( cntr_param.speed>0 && cntr_param.start==0)
483  rescale_object( wname, 0.3, SURF::Vector_3<double>(0.5, 0.75, 0.5));
484  }
485  else if ( cntr_param.test == "vortex") {
486  spd = &vortex_spd;
487  if ( cntr_param.speed>0 && cntr_param.start==0)
488  rescale_object( wname, 0.15, SURF::Vector_3<double>(0.5, 0.75, 0.5));
489  }
490  else if ( cntr_param.test == "leveque") {
491  spd = &leveque_spd;
492  if ( cntr_param.speed>0 && cntr_param.start==0)
493  rescale_object( wname, 0.15, SURF::Vector_3<double>(0.35, 0.35, 0.35));
494  }
495  else {
496  COM_assertion_msg(false, "Unkonwn test case");
497  }
498 
499  if ( cntr_param.start==0) {
500  double zero=0.;
501  // Perform one step of smoothing to detect features.
502  COM_call_function( PROP_propagate, &pmesh, &vvels, &zero, &disps, &zero);
503  }
504 
505  double vol = output_solution( wname, cntr_param.start==0?"00000":NULL);
506 
508  COM_set_profiling_barrier( PROP_propagate, MPI_COMM_WORLD);
509 
510  std::ofstream areas, vols;
511  areas.precision(10); vols.precision(10);
512  cout.precision(10);
513 
514  if ( rank==0) areas.open((wname+"_areas.m").c_str());
515  if ( rank==0) vols.open((wname+"_vols.m").c_str());
516 
517  double a=compute_area( wname), v=compute_volume(wname);
518 
519  if ( rank==0) {
520  std::cout << "Area is " << a << " and total volume is " << v << std::endl;
521  areas << wname << "_as = [0 " << a << ";..." << std::endl;
522  vols << wname << "_vs = [0 " << v << ";..." << std::endl;
523  }
524 
525  char fname[40];
526  std::sprintf(fname, "timedata_%03d.txt", rank);
527 
528 #ifdef MESH_ADAPT
529  AdaptCOMWindow acw( wname.c_str());
530 
531  acw.set_option( "adapt_iter", cntr_param.adapt_iter);
532  acw.set_option( "redist_iter", cntr_param.do_redist);
533 
534  acw.set_option( "do_collapse", cntr_param.do_collapse);
535  acw.set_option( "do_flip", cntr_param.do_flip);
536  acw.set_option( "do_split", cntr_param.do_split);
537  acw.set_option( "refine", cntr_param.refine);
538 
539  if ( !cntr_param.fangle.empty())
540  acw.set_option( "fangle_weak",
541  std::atof(cntr_param.fangle.c_str())/180*M_PI);
542  if ( !cntr_param.sangle.empty())
543  acw.set_option( "fangle_strong",
544  std::atof(cntr_param.sangle.c_str())/180*M_PI);
545  if ( !cntr_param.eigthres.empty())
546  acw.set_option( "eig_thres", std::atof(cntr_param.eigthres.c_str()));
547 
548  if ( cntr_param.collapse_ratio>0)
549  acw.set_option( "collapse_ratio", cntr_param.collapse_ratio);
550  if ( cntr_param.split_angle > 0)
551  acw.set_option( "split_angle", cntr_param.split_angle);
552 #endif
553 
554  double t=cntr_param.start*cntr_param.timestep;
555  for ( int i=1+cntr_param.start; i<=cntr_param.steps; ++i, t+=cntr_param.timestep) {
556  if ( rank==0) cout << "Step " << i << endl;
557 
558  double local_t = 0, t_rem = cntr_param.timestep, t_elapsed;
559  while ( t_rem > 0) {
560  if ( use_quadpoints) {
561  // Initialize the velocity at quadrature points.
562  PointPropagate::propagate_faces( wname, *spd, t,
563  cntr_param.timestep, "qvels");
564  }
565  else {
566  // Initialize the velocity at vertices.
567  PointPropagate::propagate_nodes( wname, *spd, t,
568  cntr_param.timestep, "vvels");
569  }
570 
571  // If speed is 0, then do smoothing only.
572  if ( cntr_param.speed==0) t_rem=0;
573 
574  COM_call_function( PROP_propagate, &pmesh, use_quadpoints?&qvels:&vvels,
575  &t_rem, &disps, &t_elapsed);
576 
577 
578  COM_call_function( BLAS_add, &nc, &disps, &nc);
579  COM_call_function( BLAS_add, &disps_total, &disps, &disps_total);
580 
581  local_t = cntr_param.timestep - (t_rem-t_elapsed);
582  if ( local_t*1.e3<cntr_param.timestep) {
583  std::cout << "Time step too small. Stopping..." << std::endl;
584  exit(-1);
585  }
586 
587 #ifdef MESH_ADAPT
588  if ( (t_elapsed<0.1*t_rem || t_elapsed==t_rem &&
589  i%cntr_param.remesh_interval==0) && cntr_param.adapt_iter>0) {
590 
591  for (int j=-1; j<cntr_param.adapt_iter; ++j) {
592  if ( j>=0) { // Note: we perform a step of smoothing first.
593  if ( acw.adapt()==0) break;
594  // If the window has changed, reinitialize Rocprop.
595  COM_call_function( PROP_init, &pmesh);
596  }
597 
598  for ( int k=0; k<cntr_param.do_redist; ++k) {
599  // Perform mesh smoothing
600  double zero=0.;
601  std::cout << "Perform anisotropic smoothing" << std::endl;
602  COM_call_function( PROP_propagate, &pmesh, &vvels,
603  &zero, &disps, &t_elapsed);
604  COM_call_function( BLAS_add, &nc, &disps, &nc);
605  COM_call_function( BLAS_add, &disps_total, &disps, &disps_total);
606  }
607  }
608  }
609 #endif
610 
611  t_rem = cntr_param.timestep - local_t;
612 
613  }
614 
615  a=compute_area( wname); v=compute_volume(wname);
616  if ( rank==0) {
617  std::cout << "Area is " << a << " and total volume is " << v << std::endl;
618  double t = i*cntr_param.timestep;
619  areas << '\t' << t << '\t';
620  areas.precision(16); areas << a << ";..." << std::endl;
621  vols << '\t' << t << '\t';
622  vols.precision(16); vols << v << ";..." << std::endl;
623  }
624 
625  if ( i%cntr_param.interval == 0) {
626  char steps[10];
627  if ( cntr_param.steps>=1.e6)
628  std::sprintf( steps, "%07d", i);
629  else if ( cntr_param.steps>=1.e5)
630  std::sprintf( steps, "%06d", i);
631  else
632  std::sprintf( steps, "%05d", i);
633 
634  output_solution( wname, steps, vol);
635 
636  COM_print_profile( fname, "Proptest");
637  }
638  }
639 
640  areas << "];" << std::endl;
641  vols << "];" << std::endl;
642 
643  COM_finalize();
644 }
645 
646 
647 
648 
649 
650 
int COMMPI_Comm_rank(MPI_Comm c)
Definition: commpi.h:162
double square(double x)
here we put it at the!beginning of the common block The point to point and collective!routines know about but MPI_TYPE_STRUCT as yet does not!MPI_STATUS_IGNORE and MPI_STATUSES_IGNORE are similar objects!Until the underlying MPI library implements the C version of these are declared as arrays of MPI_STATUS_SIZE!The types and are OPTIONAL!Their values are zero if they are not available Note that!using these reduces the portability of MPI_IO INTEGER MPI_BOTTOM INTEGER MPI_DOUBLE_PRECISION INTEGER MPI_LOGICAL INTEGER MPI_2REAL INTEGER MPI_2DOUBLE_COMPLEX INTEGER MPI_LB INTEGER MPI_WTIME_IS_GLOBAL INTEGER MPI_COMM_WORLD
static void propagate_faces(const std::string &wname, const Speed &s, double t, double dt, const std::string &vname)
j indices k indices k
Definition: Indexing.h:6
void int int REAL REAL * y
Definition: read.cpp:74
double compute_volume(const string &wname)
Definition: advectest.C:414
#define COM_assertion_msg(EX, msg)
void COM_set_size(const char *wa_str, int pane_id, int size, int ng=0)
Set sizes of for a specific attribute.
Definition: roccom_c++.h:136
This file contains the prototypes for Roccom API.
C/C++ Data types.
Definition: roccom_basic.h:129
void rescale_object(std::string &wname, double alpha, const SURF::Vector_3< double > &origin)
Definition: advectest.C:300
int COM_get_attribute_handle(const char *waname)
Definition: roccom_c++.h:412
void init_parameters(const Control_parameter &cntr_param)
Definition: advectest.C:217
Definition: speeds.h:32
float atof(const char *const str)
Read a float number from a C-string.
Definition: CImg.h:4905
*********************************************************************Illinois Open Source License ****University of Illinois NCSA **Open Source License University of Illinois All rights reserved ****Developed free of to any person **obtaining a copy of this software and associated documentation to deal with the Software without including without limitation the rights to and or **sell copies of the and to permit persons to whom the **Software is furnished to do subject to the following this list of conditions and the following disclaimers ****Redistributions in binary form must reproduce the above **copyright this list of conditions and the following **disclaimers in the documentation and or other materials **provided with the distribution ****Neither the names of the Center for Simulation of Advanced the University of nor the names of its **contributors may be used to endorse or promote products derived **from this Software without specific prior written permission ****THE SOFTWARE IS PROVIDED AS WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **ARISING OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE v
Definition: roccomf90.h:20
std::string read_in_mesh(const char *fname)
Definition: advectest.C:283
double collapse_ratio
Definition: advectest.C:110
string normaldif
Definition: advectest.C:84
void COM_finalize()
Definition: roccom_c++.h:59
static void mul_scalar(const Attribute *x, const void *y, Attribute *z, int swap=0)
Operation wrapper for multiplication with y as a scalar pointer.
Definition: op3args.C:347
void int int int REAL REAL REAL * z
Definition: write.cpp:76
double compute_area(const string &wname)
Definition: advectest.C:395
void COM_print_profile(const char *fname, const char *header)
Definition: roccom_c++.h:557
Definition: Rocout.h:81
blockLoc i
Definition: read.cpp:79
void test(void)
Definition: flotsam.C:99
void int int REAL * x
Definition: read.cpp:74
void COM_set_profiling_barrier(int hdl, MPI_Comm comm)
Definition: roccom_c++.h:554
void COM_window_init_done(const char *w_str, int pane_changed=true)
Definition: roccom_c++.h:102
This file contains a set of routines for error assertion.
string wavefrontal
Definition: advectest.C:83
void print_usage(int argc, char *argv[])
Definition: advectest.C:68
void COM_set_profiling(int i)
Definition: roccom_c++.h:550
void COM_call_function(const int wf, int argc,...)
Definition: roccom_c.C:48
void load_modules()
Definition: advectest.C:58
int main(int argc, char *argv[])
Definition: blastest.C:94
static void propagate_nodes(const std::string &wname, const Speed &s, double t, double dt, const std::string &vname)
j indices j
Definition: Indexing.h:6
void COM_init(int *argc, char ***argv)
Definition: roccom_c++.h:57
NT abs(const NT &x)
Definition: number_utils.h:130
void COM_new_attribute(const char *wa_str, const char loc, const int type, int ncomp, const char *unit)
Registering an attribute type.
Definition: roccom_c++.h:118
static int rank
Definition: advectest.C:66
COM_END_NAME_SPACE COM::Roccom_base * COM_get_roccom()
Definition: Roccom_base.h:537
double split_angle
Definition: advectest.C:111
int read_winmesh(const char *fname, const std::string &wname, bool del=true)
Definition: IM_Reader.h:58
int COMMPI_Initialized()
Definition: commpi.h:168
double output_solution(const string &wname, const char *timelevel, double ref=0.)
Definition: advectest.C:357
void init_attributes(const string &wname, const Control_parameter &cntr_param)
Definition: advectest.C:323
#define COM_LOAD_MODULE_STATIC_DYNAMIC(moduleName, windowString)
Definition: roccom_basic.h:111
void COM_resize_array(const char *wa_str, int pane_id=0, void **addr=NULL, int strd=-1, int cap=0)
Resize an attribute on a specific pane and return the address by setting addr.
Definition: roccom_c++.h:200
int COM_get_function_handle(const char *wfname)
Definition: roccom_c++.h:428
void read_control_file(const char *fname, Control_parameter &cp)
Definition: advectest.C:115
static void add_scalar(const Attribute *x, const void *y, Attribute *z, int swap=0)
Operation wrapper for addition with y as a scalar pointer.
Definition: op3args.C:299
#define COM_EXTERN_MODULE(moduleName)
Definition: roccom_basic.h:116