MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
twoTanksMacroscopicProblem.C
Go to the documentation of this file.
1 
41 #include <stdio.h>
42 #include <iostream>
43 #include "modena.h"
44 
45 using namespace std;
46 
47 int
48 main(int argc, char *argv[])
49 {
50  const double D = 0.01;
51 
52  double p0 = 3e5;
53  double p1 = 10000;
54  double V0 = 0.1;
55  double V1 = 1;
56  double T = 300;
57 
58  double t = 0.0;
59  const double deltat = 1e-3;
60  const double tend = 5.5;
61 
62  double m0 = p0*V0/287.1/T;
63  double m1 = p1*V1/287.1/T;
64 
65  double rho0 = m0/V0;
66  double rho1 = m1/V1;
67 
68 
69  // Instantiate a model
70  modena_model_t *model = modena_model_new("flowRate");
71 
72  if(modena_error_occurred())
73  {
74  return modena_error();
75  }
76 
77  // Allocate memory and fetch arg positions
78  modena_inputs_t *inputs = modena_inputs_new(model);
79  modena_outputs_t *outputs = modena_outputs_new(model);
80 
81  cout << "inputs:" << endl;
82  const char** iNames = modena_model_inputs_names(model);
83  for(int i=0; i<modena_model_inputs_size(model); i++)
84  {
85  cout << iNames[i] << endl;
86  }
87 
88  cout << "outputs:" << endl;
89  const char** oNames = modena_model_outputs_names(model);
90  for(int i=0; i<modena_model_outputs_size(model); i++)
91  {
92  cout << oNames[i] << endl;
93  }
94 
95  cout << "parameters:" << endl;
96  const char** pNames = modena_model_parameters_names(model);
97  for(int i=0; i<modena_model_parameters_size(model); i++)
98  {
99  cout << pNames[i] << endl;
100  }
101 
102  size_t Dpos = modena_model_inputs_argPos(model, "D");
103  size_t rho0Pos = modena_model_inputs_argPos(model, "rho0");
104  size_t p0Pos = modena_model_inputs_argPos(model, "p0");
105  size_t p1Byp0Pos = modena_model_inputs_argPos(model, "p1Byp0");
106 
108 
109  while(t + deltat < tend + 1e-10)
110  {
111  t += deltat;
112 
113  if(p0 > p1)
114  {
115  // Set input vector
116  modena_inputs_set(inputs, Dpos, D);
117  modena_inputs_set(inputs, rho0Pos, rho0);
118  modena_inputs_set(inputs, p0Pos, p0);
119  modena_inputs_set(inputs, p1Byp0Pos, p1/p0);
120 
121  // Call the model
122  int ret = modena_model_call(model, inputs, outputs);
123 
124  // Terminate, if requested
125  if(modena_error_occurred())
126  {
127  modena_inputs_destroy(inputs);
128  modena_outputs_destroy(outputs);
129  modena_model_destroy(model);
130 
131  return modena_error();
132  }
133 
134  // Fetch result
135  double mdot = modena_outputs_get(outputs, 0);
136 
137  m0 -= mdot*deltat;
138  m1 += mdot*deltat;
139  }
140  else
141  {
142  // Set input vector
143  modena_inputs_set(inputs, Dpos, D);
144  modena_inputs_set(inputs, rho0Pos, rho1);
145  modena_inputs_set(inputs, p0Pos, p1);
146  modena_inputs_set(inputs, p1Byp0Pos, p0/p1);
147 
148  // Call the model
149  int ret = modena_model_call(model, inputs, outputs);
150 
151  // Terminate, if requested
152  if(modena_error_occurred())
153  {
154  modena_inputs_destroy(inputs);
155  modena_outputs_destroy(outputs);
156  modena_model_destroy(model);
157 
158  return modena_error();
159  }
160 
161  // Fetch result
162  double mdot = modena_outputs_get(outputs, 0);
163 
164  m0 += mdot*deltat;
165  m1 -= mdot*deltat;
166  }
167 
168  rho0 = m0/V0;
169  rho1 = m1/V1;
170  p0 = m0/V0*287.1*T;
171  p1 = m1/V1*287.1*T;
172 
173  cout << "t = " << t << " rho0 = " << rho0 << " p0 = " << p0 << " p1 = " << p1 << endl;
174  }
175 
176  modena_inputs_destroy(inputs);
177  modena_outputs_destroy(outputs);
178  modena_model_destroy(model);
179 
180  return 0;
181 }
char ** modena_model_parameters_names(const modena_model_t *self)
Function returning the names of the parameters.
Definition: model.c:442
char ** modena_model_inputs_names(const modena_model_t *self)
Function returning the names of the inputs.
Definition: model.c:432
size_t modena_model_inputs_argPos(const modena_model_t *self, const char *name)
Function determining position of an argument in the input vector.
Definition: model.c:366
size_t modena_model_inputs_size(const modena_model_t *self)
Function returning the size of the input vector.
Definition: model.c:447
size_t modena_model_outputs_size(const modena_model_t *self)
Function returning the size of the output vector.
Definition: model.c:452
modena_model_t * modena_model_new(const char *modelId)
Function fetching a surrogate model from MongoDB.
Definition: model.c:281
int modena_model_call(modena_model_t *self, modena_inputs_t *inputs, modena_outputs_t *outputs)
Function calling the surrogate model and checking for errors.
Definition: model.c:553
int main(int argc, char *argv[])
Reads parameters. Creates struts and walls. Saves foam morphology to a file.
Definition: foams.cc:23
stores a surrogate model
Definition: model.h:94
void modena_model_destroy(modena_model_t *self)
Function deallocating the memory allocated for the surrogate model.
Definition: model.c:665
Definition: model.f90:7
size_t modena_model_parameters_size(const modena_model_t *self)
Function returning the size of the parameter vector.
Definition: model.c:457
double tend
end time, s
void modena_model_argPos_check(const modena_model_t *self)
Function checking that the user has queried all input positions.
Definition: model.c:408
char ** modena_model_outputs_names(const modena_model_t *self)
Function returning the names of the outputs.
Definition: model.c:437