MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
fullerEtAlDiffusionTest.C
1 /*
2 
3  ooo ooooo oooooooooo. ooooo ooo
4  `88. .888' `888' `Y8b `888b. `8'
5  888b d'888 .ooooo. 888 888 .ooooo. 8 `88b. 8 .oooo.
6  8 Y88. .P 888 d88' `88b 888 888 d88' `88b 8 `88b. 8 `P )88b
7  8 `888' 888 888 888 888 888 888ooo888 8 `88b.8 .oP"888
8  8 Y 888 888 888 888 d88' 888 .o 8 `888 d8( 888
9  o8o o888o `Y8bod8P' o888bood8P' `Y8bod8P' o8o `8 `Y888""8o
10 
11 Copyright
12  2014-2016 MoDeNa Consortium, All rights reserved.
13 
14 License
15  This file is part of Modena.
16 
17  Modena is free software; you can redistribute it and/or modify it under
18  the terms of the GNU General Public License as published by the Free
19  Software Foundation, either version 3 of the License, or (at your option)
20  any later version.
21 
22  Modena is distributed in the hope that it will be useful, but WITHOUT ANY
23  WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
24  FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
25  details.
26 
27  You should have received a copy of the GNU General Public License along
28  with Modena. If not, see <http://www.gnu.org/licenses/>.
29 
30 Description
31  Solving the two tank problem the MoDeNa way.
32 
33  A prototypical macros-scopic code embeds a micro-scale model (flowRate)
34  through the MoDeNa interface library.
35 
36 Authors
37  Henrik Rusche
38 
39 Contributors
40 */
41 
42 #include <stdio.h>
43 #include <iostream>
44 #include "modena.h"
45 
46 using namespace std;
47 
48 int
49 main(int argc, char *argv[])
50 {
51  {
52  // Instantiate index set
53  modena_index_set_t *indexSet = modena_index_set_new("species");
54 
55  // Find index
56  cout << "index[N2] = "
57  << modena_index_set_get_index(indexSet, "N2")
58  << endl;
59 
60  // Iterate index set and print names
61  size_t start = modena_index_set_iterator_start(indexSet);
62  size_t end = modena_index_set_iterator_end(indexSet);
63  for(int i = start; i < end; i++)
64  {
65  cout << "index[" << i << "] = "
66  << modena_index_set_get_name(indexSet, i)
67  << endl;
68  }
69 
70  modena_index_set_destroy(indexSet);
71  }
72 
73  {
74  // Instantiate function
75  modena_function_t *function =
76  modena_function_new("fullerEtAlDiffusion");
77 
78  // Get index set from function
79  modena_index_set_t *indexSet =
80  modena_function_get_index_set(function, "A");
81 
82  // Find index
83  cout << "index[N2] = "
84  << modena_index_set_get_index(indexSet, "N2")
85  << endl;
86 
87  // Iterate index set and print names
88  size_t start = modena_index_set_iterator_start(indexSet);
89  size_t end = modena_index_set_iterator_end(indexSet);
90  for(int i = start; i < end; i++)
91  {
92  cout << "index[" << i << "] = "
93  << modena_index_set_get_name(indexSet, i)
94  << endl;
95  }
96 
97  modena_index_set_destroy(indexSet);
98  // Not working - BUG!
99  //modena_function_destroy(function);
100  }
101 
102  // Instantiate model
104  (
105  "fullerEtAlDiffusion[A=H2O,B=N2]"
106  );
107 
108  if(modena_error_occurred())
109  {
110  return modena_error();
111  }
112 
113  // Allocate memory and fetch arg positions
114  modena_inputs_t *inputs = modena_inputs_new(model);
115  modena_outputs_t *outputs = modena_outputs_new(model);
116 
117  size_t ppos = modena_model_inputs_argPos(model, "p");
118  size_t TPos = modena_model_inputs_argPos(model, "T");
119 
120  size_t DPos = modena_model_outputs_argPos(model, "D[A]");
121 
123 
124  const double p = 1e5;
125  const double T = 290;
126 
127  // Set input vector
128  modena_inputs_set(inputs, ppos, p);
129  modena_inputs_set(inputs, TPos, T);
130 
131  // Call the model
132  int ret = modena_model_call(model, inputs, outputs);
133 
134  // Terminate, if requested
135  if(modena_error_occurred())
136  {
137  modena_inputs_destroy(inputs);
138  modena_outputs_destroy(outputs);
139  modena_model_destroy(model);
140 
141  return modena_error();
142  }
143 
144  // Fetch result
145  double D = modena_outputs_get(outputs, 0);
146 
147  cout << "p = " << p << " T = " << T << " D = " << D << endl;
148 
149  modena_inputs_destroy(inputs);
150  modena_outputs_destroy(outputs);
151  modena_model_destroy(model);
152 
153  return 0;
154 }
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
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_outputs_argPos(const modena_model_t *self, const char *name)
Function determining position of a result in the output vector.
Definition: model.c:392
void modena_model_argPos_check(const modena_model_t *self)
Function checking that the user has queried all input positions.
Definition: model.c:408