xref: /phasta/phSolver/incompressible/usr.c (revision c07b23fc5d98026b0b2c751adbe455d84e9bf2d1)
1 /*===========================================================================
2  *
3  * "usr.c":  user's function
4  *
5  *===========================================================================
6  */
7 #include <stdio.h>
8 #include <stdlib.h>
9 #include "les.h"
10 #include "usr.h"
11 //mr change
12 // #include "common_c.h" //not needed here any more because added in new_interface.c
13 //mr change end
14 #include "common_c.h"
15 #include "phastaIO.h"
16 #include "phIO.h"
17 #include "new_interface.h"
18 #include <FCMangle.h>
19 
20 extern char phasta_iotype[80];
21 extern int field_flag; //new_interface.c //TODO: figure out where these
22 extern int f_descriptor; //new_interface.c // should be
23 
24 /*===========================================================================
25  *
26  * "usrNew":  Put all the values in usrHd
27  *
28  * From FORTRAN
29  *
30  *	integer		usr(100)
31  *	dimension	aperm(numnp,nperm)
32  *	...
33  *	call usrnew( usr, aperm, ..., numnp, ...)
34  *
35  *
36  *===========================================================================
37  */
38 #include "mpi.h"
39 static int lmNum = 0;
40 static LesHd lesArray[8];
41 void   usrNew(	UsrHd	  usrHd,
42                         int*      eqnType,
43                         double*	  aperm,
44                         double*	  atemp,
45                         double*   resf,
46                         double*   solinc,
47                         double*   flowDiag,
48                         double*   sclrDiag,
49                         double*   lesP,
50                         double*   lesQ,
51                         int*      iBC,
52                         double*   BC,
53                         int*      iper,
54                         int*      ilwork,
55                         int*      numpe,
56                         int*      nNodes,
57                         int*      nenl,
58                         int*	  nPermDims,
59                         int*	  nTmpDims,
60                         int*	  rowp,
61                         int*	  colm,
62                         double*   lhsK,
63                         double*   lhsP,
64                         double*   lhsS,
65                         int*      nnz_tot,
66                         double*   CGsol
67     )
68 {
69     char*	funcName = "usrNew" ;	/* function name		*/
70 
71 /*---------------------------------------------------------------------------
72  * Stick the parameters
73  *---------------------------------------------------------------------------
74  */
75     usrHd->eqnType      = *eqnType ;
76     usrHd->aperm	= aperm ;
77     usrHd->atemp	= atemp ;
78     usrHd->resf         = resf ;
79     usrHd->solinc       = solinc ;
80     usrHd->flowDiag     = flowDiag ;
81     usrHd->sclrDiag     = sclrDiag ;
82     usrHd->lesP         = lesP ;
83     usrHd->lesQ         = lesQ ;
84     usrHd->iBC          = iBC  ;
85     usrHd->BC           = BC   ;
86     usrHd->iper         = iper ;
87     usrHd->ilwork       = ilwork ;
88     usrHd->numpe        = *numpe ;
89     usrHd->nNodes	= *nNodes ;
90     usrHd->nenl         = *nenl ;
91     usrHd->nPermDims	= *nPermDims ;
92     usrHd->nTmpDims	= *nTmpDims ;
93     usrHd->rowp	        = rowp ;
94     usrHd->colm	        = colm ;
95     usrHd->lhsK	        = lhsK ;
96     usrHd->lhsP	        = lhsP ;
97     usrHd->lhsS         = lhsS ;
98     usrHd->nnz_tot      = nnz_tot ;
99     usrHd->CGsol        = CGsol;
100 } /* end of usrNew() */
101 
102 /*===========================================================================
103  *
104  * "usrPointer":  Get the pointer
105  *
106  *===========================================================================
107  */
108 Real*
109 usrPointer(	UsrHd	usrHd,
110             Integer	id,
111             Integer	offset,
112             Integer	nDims )
113 {
114     char*	funcName = "usrPointer";/* function name		*/
115     Real*	pnt ;			/* pointer			*/
116 
117 /*---------------------------------------------------------------------------
118  * Get the head of the memory
119  *---------------------------------------------------------------------------
120  */
121     if ( id == LES_RES_PNT ) {
122 
123         pnt	= usrHd->resf ;
124         id	= 0 ;
125 
126     } else if ( id == LES_SOL_PNT ) {
127 
128         pnt	= usrHd->solinc ;
129         id	= 0 ;
130 
131     } else if ( id < 0 ) {
132 
133         pnt	= usrHd->aperm ;
134         id	= id + usrHd->nPermDims ;
135 
136     } else {
137 
138         pnt	= usrHd->atemp ;
139         id	= id ;
140 
141     }
142 /*---------------------------------------------------------------------------
143  * Get the offset
144  *---------------------------------------------------------------------------
145  */
146     pnt		= pnt + (id + offset) * usrHd->nNodes ;
147 
148 /*---------------------------------------------------------------------------
149  * Return the pointer
150  *---------------------------------------------------------------------------
151  */
152     return( pnt ) ;
153 
154 } /* end of usrPointer() */
155 
156 #define myflesnew FortranCInterface_GLOBAL_(myflesnew,MYFLESNEW)
157 #define myflessolve FortranCInterface_GLOBAL_(myflessolve,MYFLESSOLVE)
158 #define savelesrestart FortranCInterface_GLOBAL_(savelesrestart,SAVELESRESTART)
159 #define readlesrestart FortranCInterface_GLOBAL_(readlesrestart,READLESRESTART)
160 #define solverlicenseserver FortranCInterface_GLOBAL_(solverlicenseserver,SOLVERLICENSESERVER)
161 
162 
163 
164 #ifdef intel
165         lesArray[ *lesId ] = lesNew( fileName, *lmport, &lmNum, *eqnType,
166                                      *nDofs, *minIters, *maxIters, *nKvecs,
167                                      *prjFlag, *nPrjs, *presPrjFlag, *nPresPrjs,presPrec,
168                                      *tol, *presTol, *verbose, stats, nPermDims,
169                                      nTmpDims );
170     return ;}
171 /* the following is a fake function that was required when we moved to
172    a C++ main on in the MS Visual Studio environment.  It fails to
173    link because it is looking for this function
174 */
175 void  _CrtDbgReport() {
176     return ;}
177 
178 double __vcos_(double fg) { fflush(stdout); printf(" vcos got called \n"); fflush(stdout);}
179 double __vlog_(double fg)  { fflush(stdout); printf(" vlog got called \n"); fflush(stdout);}
180 
181 
182 #endif /* we are in unix land... whew.  secretly we have equivalenced fileName and  */
183 
184 /* #ifdef LINUX*/
185 /* void flush_(int* junk ){ return; }*/
186 /* #endif*/
187 void    myflesnew(	     Integer*	lesId,
188                          Integer*	lmport,
189                          Integer*	eqnType,
190                          Integer*	nDofs,
191                          Integer*	minIters,
192                          Integer*	maxIters,
193                          Integer*	nKvecs,
194                          Integer*	prjFlag,
195                          Integer*	nPrjs,
196                          Integer*	presPrjFlag,
197                          Integer*	nPresPrjs,
198                          Real*	    tol,
199                          Real*     	presTol,
200                          Integer*	verbose,
201                          Real*     	stats,
202                          Integer*	nPermDims,
203                          Integer*	nTmpDims,
204                          char*      lmhost          ) {
205     int procId;
206 #ifdef AMG
207     int presPrec=1;
208 #else
209     int presPrec=0;
210 #endif
211     MPI_Comm_rank( MPI_COMM_WORLD, &procId ) ;
212     if(lmNum==0){
213         if(procId==0){
214             lesArray[ *lesId ] = lesNew( lmhost, *lmport, &lmNum, *eqnType,
215                                          *nDofs, *minIters, *maxIters, *nKvecs,
216                                          *prjFlag, *nPrjs, *presPrjFlag, *nPresPrjs,presPrec,
217                                          *tol, *presTol, *verbose, stats, nPermDims,
218                                          nTmpDims );
219             MPI_Bcast( &lmNum, 1, MPI_INT, 0, MPI_COMM_WORLD ) ;
220         } else {
221             MPI_Bcast( &lmNum, 1, MPI_INT, 0, MPI_COMM_WORLD ) ;
222             lesArray[ *lesId ] = lesNew( lmhost, *lmport, &lmNum, *eqnType,
223                                          *nDofs, *minIters, *maxIters, *nKvecs,
224                                          *prjFlag, *nPrjs, *presPrjFlag, *nPresPrjs,presPrec,
225                                          *tol, *presTol, *verbose, stats, nPermDims,
226                                          nTmpDims );
227         }
228     } else {
229         lesArray[ *lesId ] = lesNew( lmhost, *lmport, &lmNum, *eqnType,
230                                      *nDofs, *minIters, *maxIters, *nKvecs,
231                                      *prjFlag, *nPrjs, *presPrjFlag, *nPresPrjs,presPrec,
232                                      *tol, *presTol, *verbose, stats, nPermDims,
233                                      nTmpDims );
234     }
235     return ;
236 }
237 
238 
239 void
240 savelesrestart( Integer* lesId,
241                  Real*    aperm,
242                  Integer* nshg,
243                  Integer* myrank,
244                  Integer* lstep,
245                  Integer* nPermDims ) {
246 
247     int nPrjs, PrjSrcId;
248     int nPresPrjs, PresPrjSrcId;
249     char filename[255];
250     int iarray[3];
251     int size, nitems;
252     double* projVec;
253     int i, j, count;
254 
255 //     sprintf( filename,"restart.%d.%d", *lstep, *myrank+1 );
256 //     openfile_( filename, "append", &fileHandle );
257 
258     nPrjs = (Integer) lesGetPar( lesArray[ *lesId ], LES_ACT_PRJS );
259     PrjSrcId = (Integer) lesGetPar( lesArray[ *lesId ], LES_PRJ_VEC_ID );
260 
261     if ( PrjSrcId < 0 ) PrjSrcId += *nPermDims;
262 
263     projVec = (double*)malloc( nPrjs * ( *nshg ) * sizeof( double ) );
264 
265     count = 0;
266     for( i = PrjSrcId; i < PrjSrcId+nPrjs; i ++ ) {
267         for( j = 0 ; j < *nshg; j++ ) {
268             projVec[ count++ ] = aperm[ (*nshg) * i + j ];
269         }
270     }
271 
272     //printf("Savelessrestart, nPrjs is %d\n",nPrjs);
273 
274     iarray[ 0 ] = *nshg;
275     iarray[ 1 ] = nPrjs;
276     nitems = 2;
277     size = (*nshg)*nPrjs;
278 
279     int name_length;
280     name_length = 18;
281     Write_Field(myrank,"a","projection vectors",&name_length, (void *)projVec,"d", nshg, &nPrjs, lstep);
282 
283     //writeheader_( &fileHandle, "projection vectors ", (void*)iarray,
284     //              &nitems, &size, "double", phasta_iotype );
285     //nitems = size;
286     //writedatablock_( &fileHandle, "projection vectors ", (void*)projVec,
287     //                 &nitems, "double", phasta_iotype );
288     free(projVec);
289 
290     /*************************************************************************/
291 
292     nPresPrjs = (Integer) lesGetPar( lesArray[ *lesId ], LES_ACT_PRES_PRJS );
293     PresPrjSrcId =(Integer)lesGetPar( lesArray[ *lesId ], LES_PRES_PRJ_VEC_ID );
294     if ( PresPrjSrcId < 0 ) PresPrjSrcId += *nPermDims;
295 
296     projVec = (double*)malloc( nPresPrjs * ( *nshg ) * sizeof( double ) );
297 
298     count = 0;
299     for( i = PresPrjSrcId; i < (PresPrjSrcId + nPresPrjs) ; i ++ ) {
300         for( j = 0 ; j < *nshg; j++ ) {
301             projVec[ count++ ] = aperm[ (*nshg) * i + j ];
302         }
303     }
304 
305     iarray[ 0 ] = *nshg;
306     iarray[ 1 ] = nPresPrjs;
307     nitems = 2;
308     size = (*nshg)*nPresPrjs;
309 
310     name_length = 27;
311     Write_Field(myrank,"a","pressure projection vectors",&name_length, projVec,"d", nshg, &nPresPrjs, lstep);
312 
313 //writeheader_( &fileHandle, "pressure projection vectors ", (void*)iarray,
314 //                  &nitems, &size, "double", phasta_iotype );
315 //    nitems = size;
316 
317 //writedatablock_( &fileHandle, "pressure projection vectors ",
318 //                     (void*)projVec, &nitems, "double", phasta_iotype );
319     free( projVec);
320 
321 //closefile_( &fileHandle, "append" );
322 }
323 
324 void
325 readlesrestart( Integer* lesId,
326                  Real*    aperm,
327                  Integer* nshg,
328                  Integer* myrank,
329                  Integer* lstep ,
330                  Integer* nPermDims ) {
331 
332     int nPrjs, PrjSrcId;
333     int nPresPrjs, PresPrjSrcId;
334     char filename[255];
335     phio_fp fileHandle = NULL;
336     int iarray[3]={-1,-1,-1};
337     int size, nitems;
338     int itwo=2;
339     int lnshg;
340     double* projVec;
341     int i,j,count;
342 
343 //MR CHANGE
344 //    sprintf( filename,"restart.%d.%d", *lstep, *myrank+1 );
345 
346 //     int nfiles=2;
347 //     int numParts=8;
348 //     int nfields=6;
349     int nfiles;
350     int nfields;
351     int numParts;
352     int nprocs;
353     int nppf;
354 
355     nfiles = outpar.nsynciofiles;
356 //    nfields = outpar.nsynciofieldsreadrestart;
357     numParts = workfc.numpe;
358     nprocs = workfc.numpe;
359 //MR CHANGE END
360 
361 //    int nppf = numParts/nfiles;
362 
363 //      MPI_Comm_size(MPI_COMM_WORLD, &numprocs);
364     // Calculate number of parts each proc deal with and where it start and end ...
365     int nppp = numParts/nprocs;        // nppp : Number of parts per proc ...
366     int startpart = *myrank * nppp +1;    // Part id from which I (myrank) start ...
367     int endpart = startpart + nppp - 1;  // Part id to which I (myrank) end ...
368 
369     sprintf(filename,"restart-dat.%d.",*lstep);
370     phio_openfile_read(filename, &nfiles, &fileHandle);
371 
372     if ( fileHandle < 0 ) return; // See phastaIO.cc for error fileHandle
373     phio_readheader(fileHandle, "projection vectors", (void*)iarray,
374                 &itwo, "integer", phasta_iotype);
375 
376     if ( iarray[0] != *nshg ) {
377         phio_closefile_read(fileHandle);
378         if(workfc.myrank==workfc.master)
379           printf("projection vectors are being initialized to zero (SAFE)\n");
380         return;
381     }
382 
383     lnshg = iarray[ 0 ] ;
384     nPrjs = iarray[ 1 ] ;
385 
386     size = (*nshg)*nPrjs;
387     projVec = (double*)malloc( size * sizeof( double ));
388 
389     phio_readdatablock(fileHandle, "projection vectors", (void*)projVec,
390                     &size, "double", phasta_iotype );
391 
392     lesSetPar( lesArray[ *lesId ], LES_ACT_PRJS, (Real) nPrjs );
393     PrjSrcId = (Integer) lesGetPar( lesArray[ *lesId ], LES_PRJ_VEC_ID );
394     if ( PrjSrcId < 0 ) PrjSrcId += *nPermDims;
395 
396     count = 0;
397     for( i = PrjSrcId; i < PrjSrcId+nPrjs; i ++ ) {
398         for( j = 0 ; j < *nshg; j++ ) {
399             aperm[ (*nshg) * i + j ] = projVec[ count++ ] ;
400         }
401     }
402 
403     free( projVec );
404 
405     /************************************************************************/
406 
407     iarray[0] = -1; iarray[1] = -1; iarray[2] = -1;
408 
409     phio_readheader(fileHandle, "pressure projection vectors", (void*)iarray,
410                  &itwo, "integer", phasta_iotype );
411 
412     lnshg = iarray[ 0 ] ;
413     nPresPrjs = iarray[ 1 ] ;
414 
415     if ( lnshg != *nshg )  {
416         phio_closefile_read(fileHandle);
417         if(workfc.myrank==workfc.master)
418           printf("pressure projection vectors are being initialized to zero (SAFE)\n");
419         return;
420     }
421 
422     size = (*nshg)*nPresPrjs;
423     projVec = (double*)malloc( size * sizeof( double ));
424 
425     phio_readdatablock(fileHandle, "pressure projection vectors", (void*)projVec,
426                     &size, "double", phasta_iotype );
427 
428     lesSetPar( lesArray[ *lesId ], LES_ACT_PRES_PRJS, (Real) nPresPrjs );
429     PresPrjSrcId=(Integer)lesGetPar( lesArray[ *lesId ], LES_PRES_PRJ_VEC_ID );
430     if ( PresPrjSrcId < 0 ) PresPrjSrcId += *nPermDims;
431 
432     count = 0;
433     for( i = PresPrjSrcId; i < PresPrjSrcId+nPresPrjs; i ++ ) {
434         for( j = 0 ; j < *nshg; j++ ) {
435             aperm[ (*nshg) * i + j ] = projVec[ count++ ] ;
436         }
437     }
438 
439     free( projVec );
440 
441     phio_closefile_read(fileHandle);
442 }
443 
444 void  myflessolve( Integer* lesId,
445                     UsrHd    usrHd){
446     lesSolve( lesArray[ *lesId ], usrHd );
447 }
448 
449 
450 int solverlicenseserver(char key[]){
451 #ifdef intel
452     strcpy(key,"C:\\cygwin\\license.dat");
453 #else
454     char* env_server_name;
455     env_server_name = getenv("LES_LICENSE_SERVER");
456     if(env_server_name) strcpy(key, env_server_name);
457     else {
458         if(workfc.myrank==workfc.master) {
459           fprintf(stderr,"environment variable LES_LICENSE_SERVER not defined \n");
460           fprintf(stderr,"using wesley as default \n");
461         }
462 //MR CHANGE
463 //        strcpy(key, "wesley.scorec.rpi.edu");
464         strcpy(key, "acusim.license.scorec.rpi.edu");
465 //MR CHANGE END
466     }
467 #endif
468     return 1;
469 }
470